Public Member Functions | |
Constructors | |
| EquationSet (const Expr &eqns, const Expr &bcs, const Expr ¶ms, const Expr ¶mValues, const Array< Expr > &fields, const Array< Expr > &fieldValues) | |
| EquationSet (const Expr &eqns, const Expr &bcs, const Array< Expr > &testFunctions, const Array< Expr > &unks, const Array< Expr > &unkLinearizationPts, const Expr &unkParams, const Expr &unkParamEvalPts, const Expr ¶ms, const Expr ¶mValues, const Array< Expr > &fixedFields, const Array< Expr > &fixedFieldValues) | |
| EquationSet (const Expr &eqns, const Expr &bcs, const Array< Expr > &vars, const Array< Expr > &varLinearizationPts, const Expr ¶ms, const Expr ¶mValues, const Array< Expr > &fixedFields, const Array< Expr > &fixedFieldValues) | |
| EquationSet (const Expr &eqns, const Expr &bcs, const Array< Expr > &vars, const Array< Expr > &varLinearizationPts, const Array< Expr > &unks, const Array< Expr > &unkLinearizationPts, const Expr ¶ms, const Expr ¶mValues, const Array< Expr > &fixedFields, const Array< Expr > &fixedFieldValues) | |
| EquationSet (const RCP< FunctionSupportResolver > &fsr, const Array< Expr > &varLinearizationPts, const Array< Expr > &unkLinearizationPts, const Expr ¶mValues, const Array< Expr > &fixedFieldValues) | |
Finding integration regions for the equation set | |
| int | numRegions () const |
| const RCP< CellFilterStub > & | region (int d) const |
| int | indexForRegion (const OrderedHandle< CellFilterStub > ®ion) const |
| bool | isBCRegion (int d) const |
| const Set< OrderedHandle < CellFilterStub > > & | regionsForTestFunc (int unreducedTestID) const |
| const Set< OrderedHandle < CellFilterStub > > & | regionsForUnkFunc (int unreducedUnkID) const |
| const Array< RegionQuadCombo > & | regionQuadCombos () const |
| const Array< RegionQuadCombo > & | bcRegionQuadCombos () const |
| bool | hasVarUnkPairs (const OrderedHandle< CellFilterStub > &domain) const |
| bool | hasBCVarUnkPairs (const OrderedHandle< CellFilterStub > &domain) const |
| const RCP< Set< OrderedPair < int, int > > > & | varUnkPairs (const OrderedHandle< CellFilterStub > &domain) const |
| const RCP< Set< OrderedPair < int, int > > > & | bcVarUnkPairs (const OrderedHandle< CellFilterStub > &domain) const |
| const Expr & | expr (const RegionQuadCombo &r) const |
| const Expr & | bcExpr (const RegionQuadCombo &r) const |
| bool | hasActiveWatchFlag () const |
| int | maxWatchFlagSetting (const std::string ¶m) const |
| const RCP < FunctionSupportResolver > & | fsr () const |
Creation of evaluation context objects | |
| EvalContext | rqcToContext (ComputationType compType, const RegionQuadCombo &r) const |
| EvalContext | bcRqcToContext (ComputationType compType, const RegionQuadCombo &r) const |
Identification of RQCs to skip for given compType | |
| bool | skipRqc (ComputationType compType, const RegionQuadCombo &r) const |
| bool | skipBCRqc (ComputationType compType, const RegionQuadCombo &r) const |
Getting information about functions | |
| int | numVarBlocks () const |
| int | numUnkBlocks () const |
| int | numUnkParams () const |
| int | numFixedParams () const |
| int | numVars (int block) const |
| int | numUnks (int block) const |
| int | numVarIDs (int block) const |
| int | numUnkIDs (int block) const |
| RCP< const CommonFuncDataStub > | varFuncData (int b, int i) const |
| RCP< const CommonFuncDataStub > | unkFuncData (int b, int i) const |
| const Expr & | unkParam (int i) const |
| const Expr & | fixedParam (int i) const |
| bool | hasVarID (int fid) const |
| bool | hasUnkID (int fid) const |
| bool | hasUnkParamID (int fid) const |
| bool | hasFixedParamID (int fid) const |
| int | blockForVarID (int varID) const |
| int | blockForUnkID (int unkID) const |
Finding the functions that appear on regions | |
| const Set< int > & | varsOnRegion (int d) const |
| const Set< int > & | unksOnRegion (int d) const |
| const Set< int > & | bcVarsOnRegion (int d) const |
| const Set< int > & | bcUnksOnRegion (int d) const |
| const Array< Set< int > > & | reducedVarsOnRegion (const OrderedHandle< CellFilterStub > &r) const |
| const Array< Set< int > > & | reducedUnksOnRegion (const OrderedHandle< CellFilterStub > &r) const |
Transforming between unreduced and reduced function IDs | |
| int | reducedVarID (int varID) const |
| int | reducedUnkID (int unkID) const |
| int | reducedUnkParamID (int unkID) const |
| int | reducedFixedParamID (int unkID) const |
| int | unreducedVarID (int block, int reducedVarID) const |
| int | unreducedUnkID (int block, int reducedUnkID) const |
| int | unreducedUnkParamID (int reducedUnkParamID) const |
| int | unreducedFixedParamID (int reducedFixedParamID) const |
Information about which calculations can be done | |
| bool | isFunctionalCalculator () const |
| bool | isSensitivityCalculator () const |
| bool | hasComputationType (ComputationType compType) const |
| const Set< ComputationType > & | computationTypes () const |
Information about which functional derivatives will be computed | |
| const DerivSet & | nonzeroFunctionalDerivs (ComputationType compType, const RegionQuadCombo &r) const |
| const DerivSet & | nonzeroBCFunctionalDerivs (ComputationType compType, const RegionQuadCombo &r) const |
Private Member Functions | |
| Expr | flattenSpectral (const Expr &input) const |
| Array< Expr > | flattenSpectral (const Array< Expr > &input) const |
| void | init (const Array< Expr > &varLinearizationPts, const Array< Expr > &unkLinearizationPts, const Expr &unkParamEvalPts, const Expr &fixedParamValues, const Array< Expr > &fixedFieldValues) |
| void | addToVarUnkPairs (const OrderedHandle< CellFilterStub > &domain, const Set< int > &vars, const Set< int > &unks, const DerivSet &nonzeros, bool isBC, int verb) |
Static Private Member Functions | |
| static Expr | toList (const Array< Expr > &e) |
Private Attributes | |
| RCP< FunctionSupportResolver > | fsr_ |
| Map< OrderedHandle < CellFilterStub >, RCP< Set < OrderedPair< int, int > > > > | varUnkPairsOnRegions_ |
| Map< OrderedHandle < CellFilterStub >, RCP< Set < OrderedPair< int, int > > > > | bcVarUnkPairsOnRegions_ |
| Array< RegionQuadCombo > | regionQuadCombos_ |
| Array< RegionQuadCombo > | bcRegionQuadCombos_ |
| Map< RegionQuadCombo, Expr > | regionQuadComboExprs_ |
| Map< RegionQuadCombo, Expr > | bcRegionQuadComboExprs_ |
| Map< ComputationType, Map < RegionQuadCombo, DerivSet > > | regionQuadComboNonzeroDerivs_ |
| Map< ComputationType, Map < RegionQuadCombo, DerivSet > > | bcRegionQuadComboNonzeroDerivs_ |
| Map< ComputationType, Map < RegionQuadCombo, EvalContext > > | rqcToContext_ |
| Map< ComputationType, Map < RegionQuadCombo, EvalContext > > | bcRqcToContext_ |
| Map< ComputationType, Set < RegionQuadCombo > > | rqcToSkip_ |
| Map< ComputationType, Set < RegionQuadCombo > > | bcRqcToSkip_ |
| Array< Expr > | unkLinearizationPts_ |
| Expr | unkParamEvalPts_ |
| Expr | fixedParamEvalPts_ |
| Set< ComputationType > | compTypes_ |
| bool | isNonlinear_ |
| bool | isVariationalProblem_ |
| bool | isFunctionalCalculator_ |
| bool | isSensitivityProblem_ |
Related Functions | |
(Note that these are not member functions.) | |
| enum | ComputationType { MatrixAndVector, VectorOnly, FunctionalOnly, FunctionalAndGradient, Sensitivities } |
| Expr | Integral (const Handle< CellFilterStub > &domain, const Expr &integrand, const WatchFlag &watch=WatchFlag()) |
| Expr | Integral (const Handle< CellFilterStub > &domain, const Expr &integrand, const Handle< QuadratureFamilyStub > &quad, const WatchFlag &watch=WatchFlag()) |
| Expr | Integral (const Handle< CellFilterStub > &domain, const Expr &integrand, const Handle< QuadratureFamilyStub > &quad, const ParametrizedCurve &curve, const WatchFlag &watch=WatchFlag()) |
Source: SundanceEquationSet.cpp
Header: SundanceEquationSet.hpp
EquationSet is an object in which the symbolic specification of a problem or functional, its BCs, its test and unknown functions, and the point about which it is to be linearized are all gathered. With this information we can compile lists of which functions are defined on which regions of the domain, which is what is required for the building of DOF maps. We can't build the DOF map here because in the Sundance core we know nothing of the mesh, so we provide accessors to the information collected by the EquationSet.
This is NOT normally a user-level object. However, EquationSet is one of the most important classes for the inner workings of Sundance, so it is critical for a developer to understand it. It is used internally in the operation of user-level classes such as LinearProblem, NonlinearProblem, and Functional.
There are several modes in which one might construct an equation set. The first is where one has written out a weak form in terms of test functions. The second is where one is taking variations of some functional.
Note that "EquationSet" is a bit of a misnomer. It was originally written to deal with setting up forward problems, but it has since been extended to encompass functionals and variations. The name persists for historical reasons; there is no particular need to change it.
Weak equations or functionals are written in terms of integrals; the regions on which integration is done must be defined somehow. Because the symbolic core knows nothing of how geometry is represented in whatever frameworks it's interacting with, the region of integration can be represented only with stub classes.
In a multivariable problem it may be useful to group variables into blocks; for example, in a segregated Navier-Stokes preconditioner the linearized equations are set up as a block system with the velocities and the pressure put into different blocks:
We use the following convention for specifying block structure: variables aggregated by Expr's listing operations are considered to be within a single block. The Teuchos Array object is then used to aggregate multiple blocks.
The EquationSet class can be used to define a functional, and one can then take variations of that functional with respect to some subset of the unknown functions appearing in the functional. We'll refer to these as the variational functions. For each variational function it is necessary to specify an evaluation point, which is simply the value about which variations are taken.
This variational capability can be used to take gradients in an optimization problem, or to derive state or adjoint equations.
Some variables may remain fixed when variations are taken. For example, in PDE-constrained optimization, the state equations are derived by taking variations of the Lagrangian with respect to the adjoint variables, holding the state and design variables constant.
Every field variable given to an equation set must also be given an evaluation point. The evaluation point is another expression, which must be of one of two types: a discrete function (subtype of DiscreteFunctionStub) or a zero expression.
It is important to understand how values of evaluation points are updated. This is NOT done by rebuilding the EquationSet object with new evaluation points. Rather, it is done by resetting the functions' internal data; because the EquationSet has stored shallow copies of the evaluation points, the EquationSet is updated automatically to follow any external changes.
Every symbolic (i.e., test or unknown) function and unknown parameter has a unique integer ID known as its function ID, or funcID for short. This ID remains associated with the function, never changing, throughout the life of the function. These IDs need not be contiguous nor ordered (they are, however, guaranteed to be unique).
In building an EquationSet, we will also create other ID numbers for each function based on the position of each function within the lists of functions given as input arguments to the equation set ctor. These IDs are contiguous and ordered, with the ordering defined by position in the input argument list. We will call these "reduced IDs." The ID intrinsic to a function is called here its "unreduced ID." EquationSet provided methods for converting between reduced and unreduced IDs.
Note that a function that appears in several equation sets might have different reduced IDs in the different equation sets, but its unreduced ID will always be the same.
Definition at line 188 of file SundanceEquationSet.hpp.
| EquationSet::EquationSet | ( | const Expr & | eqns, |
| const Expr & | bcs, | ||
| const Expr & | params, | ||
| const Expr & | paramValues, | ||
| const Array< Expr > & | fields, | ||
| const Array< Expr > & | fieldValues | ||
| ) |
Set up a functional to be integrated, where all field variables are fixed to specified values. This ctor should be used when setting up a functional for evaluation without differentiation.
| eqns | The expression defining which integrals are to be done to evaluate the functional |
| bcs | The expression defining any BC-like terms that strongly replace the ordinary equations on certain subdomains. If no BC terms are appropriate for a problem, simply enter an empty expression for this argument. |
| params | Any unknown parameters appearing in the functional. Multiple parameters should be entered as a List expression. If no parameters are present, enter an empty expression. |
| paramValues | Values of the parameters en |
| fields | The field variables (i.e., unknown functions) appearing in the functional. |
| fieldValues | Evaluation points for the variables entered in the fields argument. |
Definition at line 65 of file SundanceEquationSet.cpp.
References bcRegionQuadComboNonzeroDerivs_, bcRqcToContext_, compTypes_, fsr_, Sundance::FunctionalOnly, init(), Sundance::Map< Key, Value, Compare >::put(), Sundance::Set< Key, Compare >::put(), regionQuadComboNonzeroDerivs_, rqcToContext_, and unkParamEvalPts_.
| EquationSet::EquationSet | ( | const Expr & | eqns, |
| const Expr & | bcs, | ||
| const Array< Expr > & | testFunctions, | ||
| const Array< Expr > & | unks, | ||
| const Array< Expr > & | unkLinearizationPts, | ||
| const Expr & | unkParams, | ||
| const Expr & | unkParamEvalPts, | ||
| const Expr & | params, | ||
| const Expr & | paramValues, | ||
| const Array< Expr > & | fixedFields, | ||
| const Array< Expr > & | fixedFieldValues | ||
| ) |
Set up equations written in weak form with test functions. This ctor should be used when setting up an ordinary forward problem.
| eqns | The expression defining the weak equations. This can be linear or nonlinear. |
| bcs | The expression defining any BC-like terms that strongly replace the ordinary equations on certain subdomains. If no BC terms are appropriate for a problem, simply enter an empty expression for this argument. |
| testFunctions | The test functions used in defining the weak problem. The evaluation points for these functions are zero, and need not be given as arguments. These should be subtypes of TestFunctionStub, or lists thereof. |
| unks | The unknown functions for which the weak equation will be solved. These should be subtypes of UnknownFunctionStub, or lists thereof. |
| unkLinearizationPts | The values of the unknown function about which the equations are linearized. |
| unkParams | The unknown parameters for which the weak equation will be solved. These should of type UnknownParameter, or a list thereof. |
| unkParamEvalPts | The values of the unknown parameters about which the equations are linearized. |
| params | Any parameters whose values are held fixed (i.e, not solved for). These should be of type UnknownParameter, or a list thereof. |
| paramValues | Values of the parameters entered in the params argument. |
| fixedFields | Any field variables whose values are held fixed. These should be subtypes of UnknownFunctionStub, or lists thereof. |
| fixedFieldValues | Values of the fixed field variables. argument. |
Definition at line 117 of file SundanceEquationSet.cpp.
References bcRegionQuadComboNonzeroDerivs_, bcRqcToContext_, compTypes_, flattenSpectral(), fsr_, init(), Sundance::MatrixAndVector, Sundance::Map< Key, Value, Compare >::put(), Sundance::Set< Key, Compare >::put(), regionQuadComboNonzeroDerivs_, rqcToContext_, Sundance::Sensitivities, Sundance::Expr::size(), unkParamEvalPts_, and Sundance::VectorOnly.
| EquationSet::EquationSet | ( | const Expr & | eqns, |
| const Expr & | bcs, | ||
| const Array< Expr > & | vars, | ||
| const Array< Expr > & | varLinearizationPts, | ||
| const Expr & | params, | ||
| const Expr & | paramValues, | ||
| const Array< Expr > & | fixedFields, | ||
| const Array< Expr > & | fixedFieldValues | ||
| ) |
Definition at line 253 of file SundanceEquationSet.cpp.
References bcRegionQuadComboNonzeroDerivs_, bcRqcToContext_, compTypes_, flattenSpectral(), fsr_, Sundance::FunctionalAndGradient, Sundance::FunctionalOnly, init(), isVariationalProblem_, Sundance::Map< Key, Value, Compare >::put(), Sundance::Set< Key, Compare >::put(), regionQuadComboNonzeroDerivs_, rqcToContext_, unkLinearizationPts_, and unkParamEvalPts_.
| EquationSet::EquationSet | ( | const Expr & | eqns, |
| const Expr & | bcs, | ||
| const Array< Expr > & | vars, | ||
| const Array< Expr > & | varLinearizationPts, | ||
| const Array< Expr > & | unks, | ||
| const Array< Expr > & | unkLinearizationPts, | ||
| const Expr & | params, | ||
| const Expr & | paramValues, | ||
| const Array< Expr > & | fixedFields, | ||
| const Array< Expr > & | fixedFieldValues | ||
| ) |
Set up calculation of first and second variations of a functional. This ctor should be used when deriving the linearized form of a variational problem.
Definition at line 191 of file SundanceEquationSet.cpp.
References bcRegionQuadComboNonzeroDerivs_, bcRqcToContext_, compTypes_, flattenSpectral(), fsr_, init(), isVariationalProblem_, Sundance::MatrixAndVector, Sundance::Map< Key, Value, Compare >::put(), Sundance::Set< Key, Compare >::put(), regionQuadComboNonzeroDerivs_, rqcToContext_, unkParamEvalPts_, and Sundance::VectorOnly.
| Sundance::EquationSet::EquationSet | ( | const RCP< FunctionSupportResolver > & | fsr, |
| const Array< Expr > & | varLinearizationPts, | ||
| const Array< Expr > & | unkLinearizationPts, | ||
| const Expr & | paramValues, | ||
| const Array< Expr > & | fixedFieldValues | ||
| ) |
| void EquationSet::addToVarUnkPairs | ( | const OrderedHandle< CellFilterStub > & | domain, |
| const Set< int > & | vars, | ||
| const Set< int > & | unks, | ||
| const DerivSet & | nonzeros, | ||
| bool | isBC, | ||
| int | verb | ||
| ) | [private] |
Definition at line 972 of file SundanceEquationSet.cpp.
References Sundance::Set< Key, Compare >::begin(), Sundance::Set< Key, Compare >::contains(), Sundance::Map< Key, Value, Compare >::containsKey(), Sundance::Set< Key, Compare >::end(), Sundance::Map< Key, Value, Compare >::get(), Sundance::Deriv::isFunctionalDeriv(), Sundance::MultipleDeriv::order(), Sundance::Map< Key, Value, Compare >::put(), and SUNDANCE_MSG2.
Referenced by init().
| const Expr& Sundance::EquationSet::bcExpr | ( | const RegionQuadCombo & | r | ) | const [inline] |
Returns the BC integrand on rqc r.
Definition at line 415 of file SundanceEquationSet.hpp.
References bcRegionQuadComboExprs_, and Sundance::Map< Key, Value, Compare >::get().
| const Array<RegionQuadCombo>& Sundance::EquationSet::bcRegionQuadCombos | ( | ) | const [inline] |
Returns the list of distinct subregion-quadrature combinations appearing in the boundary conditions
Definition at line 385 of file SundanceEquationSet.hpp.
References bcRegionQuadCombos_.
Referenced by hasActiveWatchFlag(), and maxWatchFlagSetting().
| EvalContext EquationSet::bcRqcToContext | ( | ComputationType | compType, |
| const RegionQuadCombo & | r | ||
| ) | const |
Map BC RQC to the context for the derivs of the given compType
Definition at line 1166 of file SundanceEquationSet.cpp.
References bcRqcToContext_, Sundance::Map< Key, Value, Compare >::containsKey(), and Sundance::Map< Key, Value, Compare >::get().
| const Set< int > & EquationSet::bcUnksOnRegion | ( | int | d | ) | const |
Returns the unknown functions that appear in BCs on the d-th region. We can use this information to tag certain columns as BC columns in the event we're doing symmetrized BC application
Definition at line 1371 of file SundanceEquationSet.cpp.
References fsr_.
| const Set< int > & EquationSet::bcVarsOnRegion | ( | int | d | ) | const |
Returns the variational functions that appear in BCs on the d-th region. We can use this information to tag certain rows as BC rows
Definition at line 1365 of file SundanceEquationSet.cpp.
References fsr_.
| const RCP< Set< OrderedPair< int, int > > > & EquationSet::bcVarUnkPairs | ( | const OrderedHandle< CellFilterStub > & | domain | ) | const |
Returns the (var, unk) pairs appearing on the given domain. This is required for determining the sparsity structure of the matrix
Definition at line 1130 of file SundanceEquationSet.cpp.
References bcVarUnkPairsOnRegions_, Sundance::Map< Key, Value, Compare >::containsKey(), and Sundance::Map< Key, Value, Compare >::get().
| int EquationSet::blockForUnkID | ( | int | unkID | ) | const |
get the block number for the unknown function having the specified unreduced funcID
Definition at line 1259 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::GrouperBase::extractWeakForm(), and Sundance::TrivialGrouper::findGroups().
| int EquationSet::blockForVarID | ( | int | varID | ) | const |
get the block number for the variational function having the specified unreduced funcID
Definition at line 1254 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::GrouperBase::extractWeakForm(), and Sundance::TrivialGrouper::findGroups().
| const Set<ComputationType>& Sundance::EquationSet::computationTypes | ( | ) | const [inline] |
Return the types of computations this object can perform
Definition at line 611 of file SundanceEquationSet.hpp.
References compTypes_.
| const Expr& Sundance::EquationSet::expr | ( | const RegionQuadCombo & | r | ) | const [inline] |
Returns the integrand on the rqc r.
Definition at line 411 of file SundanceEquationSet.hpp.
References Sundance::Map< Key, Value, Compare >::get(), and regionQuadComboExprs_.
| const Expr& Sundance::EquationSet::fixedParam | ( | int | i | ) | const |
Returns the i-th unknown parameter
| Expr EquationSet::flattenSpectral | ( | const Expr & | input | ) | const [private] |
Flatten a spectral expression into a list of its coefficients
Definition at line 1095 of file SundanceEquationSet.cpp.
References Sundance::Expr::flatten(), Playa::Handle< PointerType >::ptr(), and Sundance::Expr::size().
Referenced by EquationSet(), and flattenSpectral().
| Array< Expr > EquationSet::flattenSpectral | ( | const Array< Expr > & | input | ) | const [private] |
Flatten a spectral expression into a list of its coefficients
Definition at line 1084 of file SundanceEquationSet.cpp.
References flattenSpectral().
| const RCP<FunctionSupportResolver>& Sundance::EquationSet::fsr | ( | ) | const [inline] |
Definition at line 426 of file SundanceEquationSet.hpp.
References fsr_.
| bool EquationSet::hasActiveWatchFlag | ( | ) | const |
Indicates whether any watch flags are active
Definition at line 1053 of file SundanceEquationSet.cpp.
References bcRegionQuadCombos(), and regionQuadCombos().
Referenced by maxWatchFlagSetting().
| bool Sundance::EquationSet::hasBCVarUnkPairs | ( | const OrderedHandle< CellFilterStub > & | domain | ) | const [inline] |
Indicates whether any BC var-unk pairs appear in the given domain
Definition at line 394 of file SundanceEquationSet.hpp.
References bcVarUnkPairsOnRegions_, and Sundance::Map< Key, Value, Compare >::containsKey().
| bool Sundance::EquationSet::hasComputationType | ( | ComputationType | compType | ) | const [inline] |
Indicate whether this equation set will do the given computation type
Definition at line 607 of file SundanceEquationSet.hpp.
References compTypes_, and Sundance::Set< Key, Compare >::contains().
| bool EquationSet::hasFixedParamID | ( | int | fid | ) | const |
Determine whether a given func ID is listed as a fixed parameter in this equation set
Definition at line 1348 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::GrouperBase::extractWeakForm().
| bool EquationSet::hasUnkID | ( | int | fid | ) | const |
Determine whether a given func ID is listed as a unk function in this equation set
Definition at line 1338 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::TrivialGrouper::findGroups().
| bool EquationSet::hasUnkParamID | ( | int | fid | ) | const |
Determine whether a given func ID is listed as a unk parameter in this equation set
Definition at line 1343 of file SundanceEquationSet.cpp.
References fsr_.
| bool EquationSet::hasVarID | ( | int | fid | ) | const |
Determine whether a given func ID is listed as a variational function in this equation set
Definition at line 1333 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::GrouperBase::extractWeakForm(), and Sundance::TrivialGrouper::findGroups().
| bool Sundance::EquationSet::hasVarUnkPairs | ( | const OrderedHandle< CellFilterStub > & | domain | ) | const [inline] |
Indicates whether any var-unk pairs appear in the given domain
Definition at line 389 of file SundanceEquationSet.hpp.
References Sundance::Map< Key, Value, Compare >::containsKey(), and varUnkPairsOnRegions_.
| int EquationSet::indexForRegion | ( | const OrderedHandle< CellFilterStub > & | region | ) | const |
Returns the index of the given region
Definition at line 1274 of file SundanceEquationSet.cpp.
References fsr_.
| void EquationSet::init | ( | const Array< Expr > & | varLinearizationPts, |
| const Array< Expr > & | unkLinearizationPts, | ||
| const Expr & | unkParamEvalPts, | ||
| const Expr & | fixedParamValues, | ||
| const Array< Expr > & | fixedFieldValues | ||
| ) | [private] |
Common initialization function called by all constructors
Definition at line 318 of file SundanceEquationSet.cpp.
References addToVarUnkPairs(), bcRegionQuadComboExprs_, bcRegionQuadComboNonzeroDerivs_, bcRegionQuadCombos_, bcRqcToContext_, bcRqcToSkip_, Sundance::Set< Key, Compare >::begin(), compTypes_, Sundance::Set< Key, Compare >::contains(), Sundance::RegionQuadCombo::domain(), Sundance::Set< Key, Compare >::elements(), Sundance::Set< Key, Compare >::end(), Sundance::SumOfIntegrals::eqnSetSetupVerb(), Sundance::Expr::flatten(), fsr_, Sundance::FunctionalAndGradient, Sundance::FunctionalOnly, Sundance::SumOfIntegrals::hasWatchedTerm(), Sundance::WatchFlag::isActive(), isVariationalProblem_, Sundance::List(), Sundance::MatrixAndVector, Playa::max(), Sundance::EvalContext::nextID(), Sundance::WatchFlag::param(), Sundance::Map< Key, Value, Compare >::put(), Sundance::Set< Key, Compare >::put(), Sundance::RegionQuadCombo::quad(), regionQuadComboExprs_, regionQuadComboNonzeroDerivs_, regionQuadCombos_, rqcToContext_, Sundance::SumOfIntegrals::rqcToExprMap(), rqcToSkip_, Sundance::Sensitivities, Sundance::EvalContext::setEvalSetupVerbosity(), Sundance::EvalContext::setSetupVerbosity(), Sundance::SymbPreprocessor::setupFunctional(), Sundance::SymbPreprocessor::setupFwdProblem(), Sundance::SymbPreprocessor::setupGradient(), Sundance::SymbPreprocessor::setupSensitivities(), Sundance::SymbPreprocessor::setupVariations(), Sundance::Set< Key, Compare >::size(), SUNDANCE_BANNER1, SUNDANCE_MSG1, SUNDANCE_MSG2, toList(), Sundance::VectorOnly, Playa::ObjectWithVerbosity::verb(), and Sundance::RegionQuadCombo::watch().
Referenced by EquationSet().
| bool EquationSet::isBCRegion | ( | int | d | ) | const |
Indicate whether the given region has an essential BC expression
Definition at line 1144 of file SundanceEquationSet.cpp.
References fsr_.
| bool Sundance::EquationSet::isFunctionalCalculator | ( | ) | const [inline] |
Definition at line 600 of file SundanceEquationSet.hpp.
References isFunctionalCalculator_.
| bool Sundance::EquationSet::isSensitivityCalculator | ( | ) | const [inline] |
Definition at line 603 of file SundanceEquationSet.hpp.
References isSensitivityProblem_.
| int EquationSet::maxWatchFlagSetting | ( | const std::string & | param | ) | const |
Finds the maximum setting of the named watch type (e.g., "fill") across all terms in the equation
Definition at line 1067 of file SundanceEquationSet.cpp.
References bcRegionQuadCombos(), hasActiveWatchFlag(), and regionQuadCombos().
| const DerivSet & EquationSet::nonzeroBCFunctionalDerivs | ( | ComputationType | compType, |
| const RegionQuadCombo & | r | ||
| ) | const |
Returns the set of nonzero functional derivatives appearing in the boundary conditions at the given subregion-quadrature combination
Definition at line 1214 of file SundanceEquationSet.cpp.
References bcRegionQuadComboNonzeroDerivs_, Sundance::Map< Key, Value, Compare >::containsKey(), and Sundance::Map< Key, Value, Compare >::get().
| const DerivSet & EquationSet::nonzeroFunctionalDerivs | ( | ComputationType | compType, |
| const RegionQuadCombo & | r | ||
| ) | const |
Returns the set of nonzero functional derivatives appearing in the equation set at the given subregion-quadrature combination
Definition at line 1204 of file SundanceEquationSet.cpp.
References Sundance::Map< Key, Value, Compare >::containsKey(), Sundance::Map< Key, Value, Compare >::get(), and regionQuadComboNonzeroDerivs_.
| int EquationSet::numFixedParams | ( | ) | const |
Returns the number of fixed parameters
Definition at line 1299 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numRegions | ( | ) | const |
Returns the number of regions on which pieces of the equation or BCs are defined.
Definition at line 1279 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numUnkBlocks | ( | ) | const |
Returns the number of unknown function blocks
Definition at line 1291 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numUnkIDs | ( | int | block | ) | const |
Returns the number of unk function IDs in this block. See the comment in FSR.hpp for an explanation of the difference between this and numVars().
Definition at line 1315 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numUnkParams | ( | ) | const |
Returns the number of unknown parameters
Definition at line 1295 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numUnks | ( | int | block | ) | const |
Returns the number of unk functions in this block
Definition at line 1307 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numVarBlocks | ( | ) | const |
Returns the number of variational function blocks
Definition at line 1287 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numVarIDs | ( | int | block | ) | const |
Returns the number of variational function IDs in this block. See the comment in FSR.hpp for an explanation of the difference between this and numVars().
Definition at line 1311 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::numVars | ( | int | block | ) | const |
Returns the number of variational functions in this block
Definition at line 1303 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::reducedFixedParamID | ( | int | unkID | ) | const |
get the reduced ID for the fixed parameter having the given funcID
Definition at line 1242 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::GrouperBase::extractWeakForm().
| int EquationSet::reducedUnkID | ( | int | unkID | ) | const |
get the reduced ID for the unknown function having the given funcID
Definition at line 1231 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::GrouperBase::extractWeakForm(), and Sundance::TrivialGrouper::findGroups().
| int EquationSet::reducedUnkParamID | ( | int | unkID | ) | const |
get the reduced ID for the unk parameter having the given funcID
Definition at line 1237 of file SundanceEquationSet.cpp.
References fsr_.
| const Array< Set< int > > & EquationSet::reducedUnksOnRegion | ( | const OrderedHandle< CellFilterStub > & | r | ) | const |
Returns the reduced unknown functions that appear explicitly on the d-th region.
Definition at line 1381 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::reducedVarID | ( | int | varID | ) | const |
get the reduced ID for the variational function having the specified unreduced funcID
Definition at line 1226 of file SundanceEquationSet.cpp.
References fsr_.
Referenced by Sundance::GrouperBase::extractWeakForm(), and Sundance::TrivialGrouper::findGroups().
| const Array< Set< int > > & EquationSet::reducedVarsOnRegion | ( | const OrderedHandle< CellFilterStub > & | r | ) | const |
Returns the reduced variational functions that appear explicitly on the d-th region
Definition at line 1376 of file SundanceEquationSet.cpp.
References fsr_.
| const RCP< CellFilterStub > & EquationSet::region | ( | int | d | ) | const |
Returns the d-th region for this equation set
Definition at line 1281 of file SundanceEquationSet.cpp.
References fsr_.
| const Array<RegionQuadCombo>& Sundance::EquationSet::regionQuadCombos | ( | ) | const [inline] |
Returns the list of distinct subregion-quadrature combinations appearing in the equation set.
Definition at line 380 of file SundanceEquationSet.hpp.
References regionQuadCombos_.
Referenced by hasActiveWatchFlag(), and maxWatchFlagSetting().
| const Set< OrderedHandle< CellFilterStub > > & EquationSet::regionsForTestFunc | ( | int | unreducedTestID | ) | const |
Return the set of regions on which the specified test func appears.
Definition at line 1264 of file SundanceEquationSet.cpp.
References fsr_.
| const Set< OrderedHandle< CellFilterStub > > & EquationSet::regionsForUnkFunc | ( | int | unreducedUnkID | ) | const |
Return the set of regions on which the specified unknown func appears
Definition at line 1269 of file SundanceEquationSet.cpp.
References fsr_.
| EvalContext EquationSet::rqcToContext | ( | ComputationType | compType, |
| const RegionQuadCombo & | r | ||
| ) | const |
Map RQC to the context for the derivs of the given compType
Definition at line 1150 of file SundanceEquationSet.cpp.
References Sundance::Map< Key, Value, Compare >::containsKey(), Sundance::Map< Key, Value, Compare >::get(), and rqcToContext_.
| bool EquationSet::skipBCRqc | ( | ComputationType | compType, |
| const RegionQuadCombo & | r | ||
| ) | const |
Map BC RQC to the context for the derivs of the given compType
Definition at line 1193 of file SundanceEquationSet.cpp.
References bcRqcToSkip_, Sundance::Map< Key, Value, Compare >::containsKey(), and Sundance::Map< Key, Value, Compare >::get().
| bool EquationSet::skipRqc | ( | ComputationType | compType, |
| const RegionQuadCombo & | r | ||
| ) | const |
Map RQC to the context for the derivs of the given compType
Definition at line 1182 of file SundanceEquationSet.cpp.
References Sundance::Map< Key, Value, Compare >::containsKey(), Sundance::Map< Key, Value, Compare >::get(), and rqcToSkip_.
| Expr EquationSet::toList | ( | const Array< Expr > & | e | ) | [static, private] |
Helper that converts an array of expr to a list expression
Definition at line 1248 of file SundanceEquationSet.cpp.
Referenced by init().
| RCP< const CommonFuncDataStub > EquationSet::unkFuncData | ( | int | b, |
| int | i | ||
| ) | const |
Returns the i-th unknown function in block b
Definition at line 1324 of file SundanceEquationSet.cpp.
References fsr_.
| const Expr & EquationSet::unkParam | ( | int | i | ) | const |
Returns the i-th unknown parameter
Definition at line 1328 of file SundanceEquationSet.cpp.
References fsr_.
| const Set< int > & EquationSet::unksOnRegion | ( | int | d | ) | const |
Returns the unknown functions that appear explicitly on the d-th region.
Definition at line 1359 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::unreducedFixedParamID | ( | int | reducedParamID | ) | const |
get the unreduced funcID for a fixed parameter as specified by a reduced ID and block index
get the unreduced funcID for a fixed parameter as specified by a reduced ID
Definition at line 1404 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::unreducedUnkID | ( | int | block, |
| int | reducedUnkID | ||
| ) | const |
get the unreduced funcID for an unknown function as specified by a reduced ID and block index
Definition at line 1392 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::unreducedUnkParamID | ( | int | reducedUnkParamID | ) | const |
get the unreduced funcID for an unknown parameter as specified by a reduced ID and block index
get the unreduced funcID for an unknown parameter as specified by a reduced ID
Definition at line 1398 of file SundanceEquationSet.cpp.
References fsr_.
| int EquationSet::unreducedVarID | ( | int | block, |
| int | reducedVarID | ||
| ) | const |
get the unreduced funcID for a variational function as specified by a reduced ID and block index
Definition at line 1387 of file SundanceEquationSet.cpp.
References fsr_.
| RCP< const CommonFuncDataStub > EquationSet::varFuncData | ( | int | b, |
| int | i | ||
| ) | const |
Returns the i-th variational function in block b
Definition at line 1319 of file SundanceEquationSet.cpp.
References fsr_.
| const Set< int > & EquationSet::varsOnRegion | ( | int | d | ) | const |
Returns the variational functions that appear explicitly on the d-th region
Definition at line 1354 of file SundanceEquationSet.cpp.
References fsr_.
| const RCP<Set<OrderedPair<int, int> > >& Sundance::EquationSet::varUnkPairs | ( | const OrderedHandle< CellFilterStub > & | domain | ) | const [inline] |
Returns the (var, unk) pairs appearing on the given domain. This is required for determining the sparsity structure of the matrix
Definition at line 400 of file SundanceEquationSet.hpp.
References Sundance::Map< Key, Value, Compare >::get(), and varUnkPairsOnRegions_.
enum ComputationType [related] |
Specifier of what sort of calculation is to be done with an equation set
Definition at line 55 of file SundanceComputationType.hpp.
| Expr Integral | ( | const Handle< CellFilterStub > & | domain, |
| const Expr & | integrand, | ||
| const WatchFlag & | watch = WatchFlag() |
||
| ) | [related] |
Referenced by main().
| Expr Integral | ( | const Handle< CellFilterStub > & | domain, |
| const Expr & | integrand, | ||
| const Handle< QuadratureFamilyStub > & | quad, | ||
| const WatchFlag & | watch = WatchFlag() |
||
| ) | [related] |
| Expr Integral | ( | const Handle< CellFilterStub > & | domain, |
| const Expr & | integrand, | ||
| const Handle< QuadratureFamilyStub > & | quad, | ||
| const ParametrizedCurve & | curve, | ||
| const WatchFlag & | watch = WatchFlag() |
||
| ) | [related] |
Definition at line 688 of file SundanceEquationSet.hpp.
Map<ComputationType, Map<RegionQuadCombo, DerivSet> > Sundance::EquationSet::bcRegionQuadComboNonzeroDerivs_ [private] |
List of the sets of nonzero functional derivatives at each regionQuadCombo
Definition at line 696 of file SundanceEquationSet.hpp.
Referenced by EquationSet(), init(), and nonzeroBCFunctionalDerivs().
Array<RegionQuadCombo> Sundance::EquationSet::bcRegionQuadCombos_ [private] |
Definition at line 682 of file SundanceEquationSet.hpp.
Referenced by bcRegionQuadCombos(), and init().
Map<ComputationType, Map<RegionQuadCombo, EvalContext> > Sundance::EquationSet::bcRqcToContext_ [private] |
List of the contexts for each BC regionQuadCombo
Definition at line 704 of file SundanceEquationSet.hpp.
Referenced by bcRqcToContext(), EquationSet(), and init().
Map<ComputationType, Set<RegionQuadCombo> > Sundance::EquationSet::bcRqcToSkip_ [private] |
List the BC RQCs that should be skipped for each computation type.
Definition at line 712 of file SundanceEquationSet.hpp.
Referenced by init(), and skipBCRqc().
Map<OrderedHandle<CellFilterStub>, RCP<Set<OrderedPair<int, int> > > > Sundance::EquationSet::bcVarUnkPairsOnRegions_ [private] |
Map from cell filter to pairs of (varID, unkID) appearing on those cells. This is needed to construct the sparsity pattern of the matrix.
Definition at line 676 of file SundanceEquationSet.hpp.
Referenced by bcVarUnkPairs(), and hasBCVarUnkPairs().
Set<ComputationType> Sundance::EquationSet::compTypes_ [private] |
Set of the computation types supported here
Definition at line 726 of file SundanceEquationSet.hpp.
Referenced by computationTypes(), EquationSet(), hasComputationType(), and init().
fixed parameter evaluation points for this equation set
Definition at line 723 of file SundanceEquationSet.hpp.
RCP<FunctionSupportResolver> Sundance::EquationSet::fsr_ [private] |
The FunctionSupportResolver deals with associating functions with subdomains
Definition at line 666 of file SundanceEquationSet.hpp.
Referenced by bcUnksOnRegion(), bcVarsOnRegion(), blockForUnkID(), blockForVarID(), EquationSet(), fsr(), hasFixedParamID(), hasUnkID(), hasUnkParamID(), hasVarID(), indexForRegion(), init(), isBCRegion(), numFixedParams(), numRegions(), numUnkBlocks(), numUnkIDs(), numUnkParams(), numUnks(), numVarBlocks(), numVarIDs(), numVars(), reducedFixedParamID(), reducedUnkID(), reducedUnkParamID(), reducedUnksOnRegion(), reducedVarID(), reducedVarsOnRegion(), region(), regionsForTestFunc(), regionsForUnkFunc(), unkFuncData(), unkParam(), unksOnRegion(), unreducedFixedParamID(), unreducedUnkID(), unreducedUnkParamID(), unreducedVarID(), varFuncData(), and varsOnRegion().
bool Sundance::EquationSet::isFunctionalCalculator_ [private] |
Flag indicating whether this equation set is a functional calculator
Definition at line 738 of file SundanceEquationSet.hpp.
Referenced by isFunctionalCalculator().
bool Sundance::EquationSet::isNonlinear_ [private] |
Flag indicating whether this equation set is nonlinear
Definition at line 730 of file SundanceEquationSet.hpp.
bool Sundance::EquationSet::isSensitivityProblem_ [private] |
Flag indicating whether this equation set is a sensitivity problem
Definition at line 742 of file SundanceEquationSet.hpp.
Referenced by isSensitivityCalculator().
bool Sundance::EquationSet::isVariationalProblem_ [private] |
Flag indicating whether this equation set is a variational problem
Definition at line 734 of file SundanceEquationSet.hpp.
Referenced by EquationSet(), and init().
Definition at line 685 of file SundanceEquationSet.hpp.
Map<ComputationType, Map<RegionQuadCombo, DerivSet> > Sundance::EquationSet::regionQuadComboNonzeroDerivs_ [private] |
List of the sets of nonzero functional derivatives at each regionQuadCombo
Definition at line 692 of file SundanceEquationSet.hpp.
Referenced by EquationSet(), init(), and nonzeroFunctionalDerivs().
Array<RegionQuadCombo> Sundance::EquationSet::regionQuadCombos_ [private] |
Definition at line 679 of file SundanceEquationSet.hpp.
Referenced by init(), and regionQuadCombos().
Map<ComputationType, Map<RegionQuadCombo, EvalContext> > Sundance::EquationSet::rqcToContext_ [private] |
List of the contexts for each regionQuadCombo
Definition at line 700 of file SundanceEquationSet.hpp.
Referenced by EquationSet(), init(), and rqcToContext().
Map<ComputationType, Set<RegionQuadCombo> > Sundance::EquationSet::rqcToSkip_ [private] |
List the RQCs that should be skipped for each computation type. This is needed in cases such as when a domain has matrix terms but no vector terms.
Definition at line 709 of file SundanceEquationSet.hpp.
Array<Expr> Sundance::EquationSet::unkLinearizationPts_ [private] |
The point in function space about which the equations are linearized
Definition at line 717 of file SundanceEquationSet.hpp.
Referenced by EquationSet().
Expr Sundance::EquationSet::unkParamEvalPts_ [private] |
unknown parameter evaluation points for this equation set
Definition at line 720 of file SundanceEquationSet.hpp.
Referenced by EquationSet().
Map<OrderedHandle<CellFilterStub>, RCP<Set<OrderedPair<int, int> > > > Sundance::EquationSet::varUnkPairsOnRegions_ [private] |
Map from cell filter to pairs of (varID, unkID) appearing on those cells. This is needed to construct the sparsity pattern of the matrix.
Definition at line 671 of file SundanceEquationSet.hpp.
Referenced by hasVarUnkPairs(), and varUnkPairs().