00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 #include "SundanceEquationSet.hpp"
00032 #include "SundanceSymbPreprocessor.hpp"
00033 #include "SundanceFunctionSupportResolver.hpp"
00034 #include "SundanceUnknownFuncElement.hpp"
00035 #include "SundanceSpectralExpr.hpp"
00036 #include "SundanceUnknownParameterElement.hpp"
00037 #include "SundanceTestFuncElement.hpp"
00038 #include "SundanceExceptions.hpp"
00039 #include "SundanceIntegral.hpp"
00040 #include "SundanceListExpr.hpp"
00041 #include "SundanceEssentialBC.hpp"
00042 #include "SundanceSumOfIntegrals.hpp"
00043 #include "SundanceSumOfBCs.hpp"
00044 #include "SundanceOut.hpp"
00045 #include "SundanceTabs.hpp"
00046
00047
00048
00049
00050
00051 using namespace Sundance;
00052 using namespace Teuchos;
00053
00054 EquationSet::EquationSet(const Expr& eqns,
00055 const Expr& bcs,
00056 const Expr& params,
00057 const Expr& paramEvalPts,
00058 const Array<Expr>& fields,
00059 const Array<Expr>& fieldValues)
00060 : fsr_(),
00061 varUnkPairsOnRegions_(),
00062 bcVarUnkPairsOnRegions_(),
00063 regionQuadCombos_(),
00064 bcRegionQuadCombos_(),
00065 regionQuadComboExprs_(),
00066 bcRegionQuadComboExprs_(),
00067 regionQuadComboNonzeroDerivs_(),
00068 bcRegionQuadComboNonzeroDerivs_(),
00069 rqcToContext_(),
00070 bcRqcToContext_(),
00071 rqcToSkip_(),
00072 bcRqcToSkip_(),
00073 unkLinearizationPts_(),
00074 unkParamEvalPts_(),
00075 fixedParamEvalPts_(paramEvalPts),
00076 compTypes_(),
00077 isNonlinear_(false),
00078 isVariationalProblem_(true),
00079 isFunctionalCalculator_(true),
00080 isSensitivityProblem_(false)
00081 {
00082 Array<Expr> unks;
00083 Array<Expr> unkEvalPt;
00084 Array<Expr> vars;
00085 Array<Expr> varEvalPt;
00086 Expr unkParams;
00087
00088 fsr_ = rcp(new FunctionSupportResolver(eqns, bcs, vars, unks, unkParams,
00089 params, fields, true));
00090
00091 compTypes_.put(FunctionalOnly);
00092
00093 rqcToContext_.put(FunctionalOnly, Map<RegionQuadCombo, EvalContext>());
00094 bcRqcToContext_.put(FunctionalOnly, Map<RegionQuadCombo, EvalContext>());
00095
00096 regionQuadComboNonzeroDerivs_.put(FunctionalOnly,
00097 Map<RegionQuadCombo, DerivSet>());
00098 bcRegionQuadComboNonzeroDerivs_.put(FunctionalOnly,
00099 Map<RegionQuadCombo, DerivSet>());
00100
00101
00102
00103 init(varEvalPt, unkEvalPt, unkParamEvalPts_, paramEvalPts, fieldValues);
00104 }
00105
00106 EquationSet::EquationSet(const Expr& eqns,
00107 const Expr& bcs,
00108 const Array<Expr>& vars,
00109 const Array<Expr>& unks,
00110 const Array<Expr>& unkLinearizationPts,
00111 const Expr& unkParams,
00112 const Expr& unkParamEvalPts,
00113 const Expr& params,
00114 const Expr& paramEvalPts,
00115 const Array<Expr>& fixedFields,
00116 const Array<Expr>& fixedFieldValues)
00117 : fsr_(),
00118 varUnkPairsOnRegions_(),
00119 bcVarUnkPairsOnRegions_(),
00120 regionQuadCombos_(),
00121 bcRegionQuadCombos_(),
00122 regionQuadComboExprs_(),
00123 bcRegionQuadComboExprs_(),
00124 regionQuadComboNonzeroDerivs_(),
00125 bcRegionQuadComboNonzeroDerivs_(),
00126 rqcToContext_(),
00127 bcRqcToContext_(),
00128 rqcToSkip_(),
00129 bcRqcToSkip_(),
00130 unkLinearizationPts_(flattenSpectral(unkLinearizationPts)),
00131 unkParamEvalPts_(unkParamEvalPts),
00132 fixedParamEvalPts_(paramEvalPts),
00133 compTypes_(),
00134 isNonlinear_(false),
00135 isVariationalProblem_(false),
00136 isFunctionalCalculator_(false),
00137 isSensitivityProblem_(unkParams.size() > 0)
00138 {
00139 compTypes_.put(MatrixAndVector);
00140 compTypes_.put(VectorOnly);
00141
00142 fsr_ = rcp(new FunctionSupportResolver(eqns, bcs, vars,
00143 unks, unkParams,
00144 params, fixedFields, false));
00145
00146 rqcToContext_.put(MatrixAndVector, Map<RegionQuadCombo, EvalContext>());
00147 bcRqcToContext_.put(MatrixAndVector, Map<RegionQuadCombo, EvalContext>());
00148
00149 rqcToContext_.put(VectorOnly, Map<RegionQuadCombo, EvalContext>());
00150 bcRqcToContext_.put(VectorOnly, Map<RegionQuadCombo, EvalContext>());
00151
00152 regionQuadComboNonzeroDerivs_.put(MatrixAndVector,
00153 Map<RegionQuadCombo, DerivSet>());
00154 bcRegionQuadComboNonzeroDerivs_.put(MatrixAndVector,
00155 Map<RegionQuadCombo, DerivSet>());
00156
00157 regionQuadComboNonzeroDerivs_.put(VectorOnly,
00158 Map<RegionQuadCombo, DerivSet>());
00159 bcRegionQuadComboNonzeroDerivs_.put(VectorOnly,
00160 Map<RegionQuadCombo, DerivSet>());
00161
00162 if (params.size() > 0)
00163 {
00164 compTypes_.put(Sensitivities);
00165 rqcToContext_.put(Sensitivities, Map<RegionQuadCombo, EvalContext>());
00166 bcRqcToContext_.put(Sensitivities, Map<RegionQuadCombo, EvalContext>());
00167 regionQuadComboNonzeroDerivs_.put(Sensitivities,
00168 Map<RegionQuadCombo, DerivSet>());
00169 bcRegionQuadComboNonzeroDerivs_.put(Sensitivities,
00170 Map<RegionQuadCombo, DerivSet>());
00171 }
00172
00173 Array<Expr> zero;
00174
00175 init(zero, flattenSpectral(unkLinearizationPts), unkParamEvalPts_,
00176 paramEvalPts, fixedFieldValues);
00177 }
00178
00179
00180 EquationSet::EquationSet(const Expr& eqns,
00181 const Expr& bcs,
00182 const Array<Expr>& vars,
00183 const Array<Expr>& varLinearizationPts,
00184 const Array<Expr>& unks,
00185 const Array<Expr>& unkLinearizationPts,
00186 const Expr& params,
00187 const Expr& paramEvalPts,
00188 const Array<Expr>& fixedFields,
00189 const Array<Expr>& fixedFieldValues)
00190 : fsr_(),
00191 varUnkPairsOnRegions_(),
00192 bcVarUnkPairsOnRegions_(),
00193 regionQuadCombos_(),
00194 bcRegionQuadCombos_(),
00195 regionQuadComboExprs_(),
00196 bcRegionQuadComboExprs_(),
00197 regionQuadComboNonzeroDerivs_(),
00198 bcRegionQuadComboNonzeroDerivs_(),
00199 rqcToContext_(),
00200 bcRqcToContext_(),
00201 rqcToSkip_(),
00202 bcRqcToSkip_(),
00203 unkLinearizationPts_(flattenSpectral(unkLinearizationPts)),
00204 unkParamEvalPts_(),
00205 fixedParamEvalPts_(paramEvalPts),
00206 compTypes_(),
00207 isNonlinear_(false),
00208 isVariationalProblem_(true),
00209 isFunctionalCalculator_(false),
00210 isSensitivityProblem_(false)
00211 {
00212 Expr unkParams;
00213 fsr_ = rcp(new FunctionSupportResolver(eqns, bcs, vars,
00214 unks, unkParams, params, fixedFields,
00215 isVariationalProblem_));
00216
00217 compTypes_.put(MatrixAndVector);
00218 compTypes_.put(VectorOnly);
00219
00220 rqcToContext_.put(MatrixAndVector, Map<RegionQuadCombo, EvalContext>());
00221 bcRqcToContext_.put(MatrixAndVector, Map<RegionQuadCombo, EvalContext>());
00222
00223 rqcToContext_.put(VectorOnly, Map<RegionQuadCombo, EvalContext>());
00224 bcRqcToContext_.put(VectorOnly, Map<RegionQuadCombo, EvalContext>());
00225
00226 regionQuadComboNonzeroDerivs_.put(MatrixAndVector,
00227 Map<RegionQuadCombo, DerivSet>());
00228 bcRegionQuadComboNonzeroDerivs_.put(MatrixAndVector,
00229 Map<RegionQuadCombo, DerivSet>());
00230
00231 regionQuadComboNonzeroDerivs_.put(VectorOnly,
00232 Map<RegionQuadCombo, DerivSet>());
00233 bcRegionQuadComboNonzeroDerivs_.put(VectorOnly,
00234 Map<RegionQuadCombo, DerivSet>());
00235
00236 init(flattenSpectral(varLinearizationPts),
00237 flattenSpectral(unkLinearizationPts),
00238 unkParamEvalPts_,
00239 paramEvalPts, fixedFieldValues);
00240 }
00241
00242 EquationSet::EquationSet(const Expr& eqns,
00243 const Expr& bcs,
00244 const Array<Expr>& vars,
00245 const Array<Expr>& varLinearizationPts,
00246 const Expr& params,
00247 const Expr& paramEvalPts,
00248 const Array<Expr>& fixedFields,
00249 const Array<Expr>& fixedFieldValues)
00250 : fsr_(),
00251 varUnkPairsOnRegions_(),
00252 bcVarUnkPairsOnRegions_(),
00253 regionQuadCombos_(),
00254 bcRegionQuadCombos_(),
00255 regionQuadComboExprs_(),
00256 bcRegionQuadComboExprs_(),
00257 regionQuadComboNonzeroDerivs_(),
00258 bcRegionQuadComboNonzeroDerivs_(),
00259 rqcToContext_(),
00260 bcRqcToContext_(),
00261 rqcToSkip_(),
00262 bcRqcToSkip_(),
00263 unkLinearizationPts_(),
00264 unkParamEvalPts_(),
00265 fixedParamEvalPts_(paramEvalPts),
00266 compTypes_(),
00267 isNonlinear_(false),
00268 isVariationalProblem_(true),
00269 isFunctionalCalculator_(true),
00270 isSensitivityProblem_(false)
00271 {
00272 compTypes_.put(FunctionalOnly);
00273 compTypes_.put(FunctionalAndGradient);
00274
00275 Expr unkParams;
00276 Array<Expr> unks;
00277 fsr_ = rcp(new FunctionSupportResolver(eqns, bcs, vars,
00278 unks, unkParams,
00279 params, fixedFields, isVariationalProblem_));
00280
00281 rqcToContext_.put(FunctionalAndGradient, Map<RegionQuadCombo, EvalContext>());
00282 bcRqcToContext_.put(FunctionalAndGradient, Map<RegionQuadCombo, EvalContext>());
00283
00284 rqcToContext_.put(FunctionalOnly, Map<RegionQuadCombo, EvalContext>());
00285 bcRqcToContext_.put(FunctionalOnly, Map<RegionQuadCombo, EvalContext>());
00286
00287 regionQuadComboNonzeroDerivs_.put(FunctionalAndGradient,
00288 Map<RegionQuadCombo, DerivSet>());
00289 bcRegionQuadComboNonzeroDerivs_.put(FunctionalAndGradient,
00290 Map<RegionQuadCombo, DerivSet>());
00291
00292 regionQuadComboNonzeroDerivs_.put(FunctionalOnly,
00293 Map<RegionQuadCombo, DerivSet>());
00294 bcRegionQuadComboNonzeroDerivs_.put(FunctionalOnly,
00295 Map<RegionQuadCombo, DerivSet>());
00296
00297
00298 init(flattenSpectral(varLinearizationPts),
00299 flattenSpectral(unkLinearizationPts_),
00300 unkParamEvalPts_,
00301 paramEvalPts, fixedFieldValues);
00302 }
00303
00304
00305
00306
00307 void EquationSet::init(
00308 const Array<Expr>& varLinearizationPts,
00309 const Array<Expr>& unkLinearizationPts,
00310 const Expr& unkParamEvalPts,
00311 const Expr& fixedParamEvalPts,
00312 const Array<Expr>& fixedFieldValues)
00313 {
00314 Tabs tab0(0);
00315 Tabs tab1;
00316 bool hasBCs = fsr_->hasBCs();
00317 const SumOfBCs* bcSum = fsr_->bcSum();
00318 const SumOfIntegrals* integralSum = fsr_->integralSum();
00319
00320 int verb = 0;
00321
00322
00323 if (integralSum->hasWatchedTerm() || (hasBCs && bcSum->hasWatchedTerm()))
00324 {
00325 int v1 = integralSum->eqnSetSetupVerb();
00326 int v2 = 0;
00327 if (hasBCs) v2 = bcSum->eqnSetSetupVerb();
00328 verb = std::max(v1, v2);
00329 }
00330 SUNDANCE_BANNER1(verb, tab0, "EquationSet setup");
00331
00332
00333
00334 const Array<Expr>& unks = fsr_->unks();
00335 const Array<Expr>& vars = fsr_->vars();
00336 const Expr& unkParams = fsr_->unkParams();
00337 const Expr& fixedParams = fsr_->fixedParams();
00338 const Array<Expr>& fixedFields = fsr_->fixedFields();
00339
00340 SUNDANCE_MSG1(verb, tab0 << "fixed params = " << fixedParams);
00341
00342 const Set<int>& varFuncSet = fsr_->varFuncSet();
00343 const Set<int>& unkFuncSet = fsr_->unkFuncSet();
00344
00345 Set<RegionQuadCombo> rqcSet;
00346 Set<RegionQuadCombo> rqcBCSet;
00347
00348
00349
00350 Array<int> contextID = tuple(EvalContext::nextID(),
00351 EvalContext::nextID(),
00352 EvalContext::nextID(),
00353 EvalContext::nextID(),
00354 EvalContext::nextID());
00355
00356
00357
00358 for (Set<ComputationType>::const_iterator
00359 i=compTypes_.begin(); i!=compTypes_.end(); i++)
00360 {
00361 rqcToSkip_[*i] = Set<RegionQuadCombo>();
00362 bcRqcToSkip_[*i] = Set<RegionQuadCombo>();
00363 }
00364
00365 SUNDANCE_MSG1(verb, tab0 << "computation types = " << compTypes_);
00366
00367
00368
00369
00370
00371
00372
00373 SUNDANCE_MSG1(verb, tab1 << "processing integral terms");
00374 for (Sundance::Map<RegionQuadCombo, Expr>::const_iterator
00375 r=integralSum->rqcToExprMap().begin();
00376 r!=integralSum->rqcToExprMap().end(); r++)
00377 {
00378 Tabs tab15;
00379 Tabs tab2;
00380 RegionQuadCombo rqc = r->first;
00381 int rqcVerb = verb;
00382 int symbVerb = 0;
00383 if (rqc.watch().isActive())
00384 {
00385 symbVerb = rqc.watch().param("symbolic preprocessing");
00386 rqcVerb=rqc.watch().param("equation set setup");
00387 }
00388 SUNDANCE_MSG1(std::max(verb,rqcVerb), tab15 << "processing RQC = " << rqc);
00389
00390
00391 rqcSet.put(rqc);
00392 Expr term = r->second;
00393 OrderedHandle<CellFilterStub> reg = rqc.domain();
00394 OrderedHandle<QuadratureFamilyStub> quad = rqc.quad();
00395
00396 regionQuadComboExprs_.put(rqc, term);
00397
00398
00399 if (compTypes_.contains(MatrixAndVector))
00400 {
00401 SUNDANCE_MSG2(rqcVerb, tab2 << "preparing matrix/vector calculation");
00402 Tabs tab3;
00403 EvalContext context(rqc, makeSet(1,2), contextID[0]);
00404 context.setSetupVerbosity(symbVerb);
00405 DerivSet nonzeros;
00406
00407 if (isVariationalProblem_)
00408 {
00409 nonzeros = SymbPreprocessor
00410 ::setupVariations(term,
00411 toList(vars),
00412 toList(varLinearizationPts),
00413 toList(unks),
00414 toList(unkLinearizationPts),
00415 unkParams,
00416 unkParamEvalPts,
00417 toList(fixedFields),
00418 toList(fixedFieldValues),
00419 fixedParams,
00420 fixedParamEvalPts,
00421 context,
00422 MatrixAndVector);
00423 }
00424 else
00425 {
00426 nonzeros = SymbPreprocessor
00427 ::setupFwdProblem(term, toList(vars),
00428 toList(unks),
00429 toList(unkLinearizationPts),
00430 unkParams,
00431 unkParamEvalPts,
00432 fixedParams,
00433 fixedParamEvalPts,
00434 toList(fixedFields),
00435 toList(fixedFieldValues),
00436 context,
00437 MatrixAndVector);
00438 }
00439 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00440 if (nonzeros.size()==0)
00441 {
00442 rqcToSkip_[MatrixAndVector].put(rqc);
00443 }
00444 else
00445 {
00446 addToVarUnkPairs(rqc.domain(), varFuncSet, unkFuncSet,
00447 nonzeros, false, rqcVerb);
00448 rqcToContext_[MatrixAndVector].put(rqc, context);
00449 regionQuadComboNonzeroDerivs_[MatrixAndVector].put(rqc,
00450 nonzeros);
00451 }
00452 }
00453
00454
00455
00456
00457 if (compTypes_.contains(VectorOnly))
00458 {
00459 SUNDANCE_MSG2(rqcVerb, tab2 << "preparing vector-only calculation");
00460 Tabs tab3;
00461 EvalContext context(rqc, makeSet(1), contextID[1]);
00462 context.setSetupVerbosity(symbVerb);
00463 DerivSet nonzeros;
00464 if (isVariationalProblem_)
00465 {
00466 nonzeros = SymbPreprocessor
00467 ::setupVariations(term, toList(vars),
00468 toList(varLinearizationPts),
00469 toList(unks),
00470 toList(unkLinearizationPts),
00471 unkParams,
00472 unkParamEvalPts,
00473 toList(fixedFields),
00474 toList(fixedFieldValues),
00475 fixedParams,
00476 fixedParamEvalPts,
00477 context,
00478 VectorOnly);
00479 }
00480 else
00481 {
00482 nonzeros = SymbPreprocessor
00483 ::setupFwdProblem(term, toList(vars),
00484 toList(unks),
00485 toList(unkLinearizationPts),
00486 unkParams,
00487 unkParamEvalPts,
00488 fixedParams,
00489 fixedParamEvalPts,
00490 toList(fixedFields),
00491 toList(fixedFieldValues),
00492 context,
00493 VectorOnly);
00494 }
00495 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00496 if (nonzeros.size()==0)
00497 {
00498 rqcToSkip_[VectorOnly].put(rqc);
00499 }
00500 else
00501 {
00502 rqcToContext_[VectorOnly].put(rqc, context);
00503 regionQuadComboNonzeroDerivs_[VectorOnly].put(rqc, nonzeros);
00504 }
00505 }
00506
00507
00508
00509 if (compTypes_.contains(Sensitivities))
00510 {
00511 SUNDANCE_MSG2(rqcVerb, tab2 << "preparing sensitivity calculation");
00512 Tabs tab3;
00513 EvalContext context(rqc, makeSet(2), contextID[4]);
00514 context.setSetupVerbosity(symbVerb);
00515 DerivSet nonzeros;
00516 nonzeros = SymbPreprocessor
00517 ::setupSensitivities(term, toList(vars),
00518 toList(unks),
00519 toList(unkLinearizationPts),
00520 unkParams,
00521 unkParamEvalPts,
00522 fixedParams,
00523 fixedParamEvalPts,
00524 toList(fixedFields),
00525 toList(fixedFieldValues),
00526 context,
00527 Sensitivities);
00528 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00529 if (nonzeros.size()==0)
00530 {
00531 rqcToSkip_[Sensitivities].put(rqc);
00532 }
00533 else
00534 {
00535 rqcToContext_[Sensitivities].put(rqc, context);
00536 regionQuadComboNonzeroDerivs_[Sensitivities].put(rqc, nonzeros);
00537 }
00538 }
00539
00540
00541
00542
00543 if (compTypes_.contains(FunctionalOnly))
00544 {
00545 SUNDANCE_MSG2(rqcVerb, tab2 << "preparing functional calculation");
00546 Tabs tab3;
00547
00548 EvalContext context(rqc, makeSet(0), contextID[2]);
00549 context.setSetupVerbosity(symbVerb);
00550 DerivSet nonzeros;
00551 Expr fields;
00552 Expr fieldValues;
00553 if (fixedFields.size() > 0)
00554 {
00555 if (vars.size() > 0)
00556 {
00557 fields = List(toList(fixedFields), toList(vars));
00558 fields = fields.flatten();
00559 }
00560 else
00561 {
00562 fields = toList(fixedFields);
00563 }
00564 }
00565 else
00566 {
00567 if (vars.size() > 0)
00568 {
00569 fields = toList(vars);
00570 }
00571 }
00572 if (fixedFieldValues.size() > 0)
00573 {
00574 if (varLinearizationPts.size() > 0)
00575 {
00576 fieldValues = List(toList(fixedFieldValues),
00577 toList(varLinearizationPts));
00578 fieldValues = fieldValues.flatten();
00579 }
00580 else
00581 {
00582 fieldValues = toList(fixedFieldValues);
00583 }
00584 }
00585 else
00586 {
00587 if (varLinearizationPts.size() > 0)
00588 {
00589 fieldValues = toList(varLinearizationPts);
00590 }
00591 }
00592
00593 nonzeros = SymbPreprocessor
00594 ::setupFunctional(term,
00595 fixedParams,
00596 fixedParamEvalPts,
00597 fields,
00598 fieldValues,
00599 context,
00600 FunctionalOnly);
00601 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00602
00603 if (nonzeros.size()==0)
00604 {
00605 rqcToSkip_[FunctionalOnly].put(rqc);
00606 }
00607 else
00608 {
00609 rqcToContext_[FunctionalOnly].put(rqc, context);
00610 regionQuadComboNonzeroDerivs_[FunctionalOnly].put(rqc, nonzeros);
00611 }
00612 }
00613
00614 if (compTypes_.contains(FunctionalAndGradient))
00615 {
00616 SUNDANCE_MSG2(rqcVerb, tab2 << "preparing functional/gradient calculation");
00617 Tabs tab3;
00618 EvalContext context(rqc, makeSet(0,1), contextID[3]);
00619 context.setSetupVerbosity(symbVerb);
00620 DerivSet nonzeros;
00621 nonzeros = SymbPreprocessor
00622 ::setupGradient(term,
00623 toList(vars), toList(varLinearizationPts),
00624 fixedParams,
00625 fixedParamEvalPts,
00626 toList(fixedFields), toList(fixedFieldValues),
00627 context,
00628 FunctionalAndGradient);
00629
00630 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00631
00632 if (nonzeros.size()==0)
00633 {
00634 rqcToSkip_[FunctionalAndGradient].put(rqc);
00635 }
00636 else
00637 {
00638 rqcToContext_[FunctionalAndGradient].put(rqc, context);
00639 regionQuadComboNonzeroDerivs_[FunctionalAndGradient].put(rqc, nonzeros);
00640 }
00641 }
00642 }
00643
00644
00645 if (hasBCs)
00646 {
00647
00648
00649 SUNDANCE_MSG1(verb, tab1 << "processing essential BC terms");
00650 for (Sundance::Map<RegionQuadCombo, Expr>::const_iterator
00651 r=bcSum->rqcToExprMap().begin();
00652 r!=bcSum->rqcToExprMap().end(); r++)
00653 {
00654 Tabs tab15;
00655 RegionQuadCombo rqc = r->first;
00656 int rqcVerb = verb;
00657 int symbVerb = 0;
00658 if (rqc.watch().isActive())
00659 {
00660 symbVerb = rqc.watch().param("symbolic preprocessing");
00661 rqcVerb=rqc.watch().param("equation set setup");
00662 }
00663 SUNDANCE_MSG1(verb, tab15 << "processing BC RQC = " << rqc);
00664
00665 rqcBCSet.put(rqc);
00666 Expr term = r->second;
00667 OrderedHandle<CellFilterStub> reg = rqc.domain();
00668 OrderedHandle<QuadratureFamilyStub> quad = rqc.quad();
00669
00670 bcRegionQuadComboExprs_.put(rqc, term);
00671
00672
00673
00674
00675 if (compTypes_.contains(MatrixAndVector))
00676 {
00677 Tabs tab3;
00678 SUNDANCE_MSG2(rqcVerb, tab3 << "preparing matrix/vector calculation");
00679 EvalContext context(rqc, makeSet(1,2), contextID[0]);
00680 context.setSetupVerbosity(symbVerb);
00681 DerivSet nonzeros;
00682
00683 if (isVariationalProblem_)
00684 {
00685 nonzeros = SymbPreprocessor
00686 ::setupVariations(term,
00687 toList(vars),
00688 toList(varLinearizationPts),
00689 toList(unks),
00690 toList(unkLinearizationPts),
00691 unkParams,
00692 unkParamEvalPts,
00693 toList(fixedFields),
00694 toList(fixedFieldValues),
00695 fixedParams,
00696 fixedParamEvalPts,
00697 context,
00698 MatrixAndVector);
00699 }
00700 else
00701 {
00702 nonzeros = SymbPreprocessor
00703 ::setupFwdProblem(term, toList(vars), toList(unks),
00704 toList(unkLinearizationPts),
00705 unkParams,
00706 unkParamEvalPts,
00707 fixedParams,
00708 fixedParamEvalPts,
00709 toList(fixedFields),
00710 toList(fixedFieldValues),
00711 context,
00712 MatrixAndVector);
00713 }
00714 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00715 if (nonzeros.size()==0)
00716 {
00717 bcRqcToSkip_[MatrixAndVector].put(rqc);
00718 }
00719 else
00720 {
00721 addToVarUnkPairs(rqc.domain(), varFuncSet, unkFuncSet,
00722 nonzeros, true, rqcVerb);
00723 bcRqcToContext_[MatrixAndVector].put(rqc, context);
00724 bcRegionQuadComboNonzeroDerivs_[MatrixAndVector].put(rqc,
00725 nonzeros);
00726 }
00727 }
00728
00729
00730
00731
00732
00733
00734 if (compTypes_.contains(VectorOnly))
00735 {
00736 Tabs tab3;
00737 SUNDANCE_MSG2(rqcVerb, tab3 << "preparing vector-only calculation");
00738 EvalContext context(rqc, makeSet(1), contextID[1]);
00739 context.setSetupVerbosity(symbVerb);
00740 DerivSet nonzeros;
00741 if (isVariationalProblem_)
00742 {
00743 nonzeros = SymbPreprocessor
00744 ::setupVariations(term, toList(vars),
00745 toList(varLinearizationPts),
00746 toList(unks),
00747 toList(unkLinearizationPts),
00748 unkParams,
00749 unkParamEvalPts,
00750 toList(fixedFields),
00751 toList(fixedFieldValues),
00752 fixedParams,
00753 fixedParamEvalPts,
00754 context,
00755 VectorOnly);
00756 }
00757 else
00758 {
00759 nonzeros = SymbPreprocessor
00760 ::setupFwdProblem(term, toList(vars), toList(unks),
00761 toList(unkLinearizationPts),
00762 unkParams,
00763 unkParamEvalPts,
00764 fixedParams,
00765 fixedParamEvalPts,
00766 toList(fixedFields),
00767 toList(fixedFieldValues),
00768 context,
00769 VectorOnly);
00770 }
00771 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00772 if (nonzeros.size()==0)
00773 {
00774 bcRqcToSkip_[VectorOnly].put(rqc);
00775 }
00776 else
00777 {
00778 bcRqcToContext_[VectorOnly].put(rqc, context);
00779 bcRegionQuadComboNonzeroDerivs_[VectorOnly].put(rqc, nonzeros);
00780 }
00781 }
00782
00783
00784
00785
00786
00787
00788
00789
00790
00791 if (compTypes_.contains(Sensitivities))
00792 {
00793 Tabs tab3;
00794 SUNDANCE_MSG2(rqcVerb, tab3 << "preparing sensitivity calculation");
00795 EvalContext context(rqc, makeSet(2), contextID[4]);
00796 context.setSetupVerbosity(symbVerb);
00797 DerivSet nonzeros;
00798 nonzeros = SymbPreprocessor
00799 ::setupSensitivities(term, toList(vars), toList(unks),
00800 toList(unkLinearizationPts),
00801 unkParams,
00802 unkParamEvalPts,
00803 fixedParams,
00804 fixedParamEvalPts,
00805 toList(fixedFields),
00806 toList(fixedFieldValues),
00807 context,
00808 Sensitivities);
00809
00810 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00811 if (nonzeros.size()==0)
00812 {
00813 bcRqcToSkip_[Sensitivities].put(rqc);
00814 }
00815 else
00816 {
00817 bcRqcToContext_[Sensitivities].put(rqc, context);
00818 bcRegionQuadComboNonzeroDerivs_[Sensitivities].put(rqc, nonzeros);
00819 }
00820 }
00821
00822
00823
00824
00825
00826
00827
00828
00829 if (compTypes_.contains(FunctionalOnly))
00830 {
00831 Tabs tab3;
00832 SUNDANCE_MSG2(rqcVerb, tab3 << "preparing functional-only calculation");
00833 EvalContext context(rqc, makeSet(0), contextID[2]);
00834 context.setSetupVerbosity(symbVerb);
00835 DerivSet nonzeros;
00836 Expr fields;
00837 Expr fieldValues;
00838 if (fixedFields.size() > 0)
00839 {
00840 if (vars.size() > 0)
00841 {
00842 fields = List(toList(fixedFields), toList(vars));
00843 fields = fields.flatten();
00844 }
00845 else
00846 {
00847 fields = toList(fixedFields);
00848 }
00849 }
00850 else
00851 {
00852 if (vars.size() > 0)
00853 {
00854 fields = toList(vars);
00855 }
00856 }
00857 if (fixedFieldValues.size() > 0)
00858 {
00859 if (varLinearizationPts.size() > 0)
00860 {
00861 fieldValues = List(toList(fixedFieldValues),
00862 toList(varLinearizationPts));
00863 fieldValues = fieldValues.flatten();
00864 }
00865 else
00866 {
00867 fieldValues = toList(fixedFieldValues);
00868 }
00869 }
00870 else
00871 {
00872 if (varLinearizationPts.size() > 0)
00873 {
00874 fieldValues = toList(varLinearizationPts);
00875 }
00876 }
00877 nonzeros = SymbPreprocessor
00878 ::setupFunctional(term,
00879 fixedParams,
00880 fixedParamEvalPts,
00881 fields, fieldValues,
00882 context,
00883 FunctionalOnly);
00884
00885 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00886 if (nonzeros.size()==0)
00887 {
00888 bcRqcToSkip_[FunctionalOnly].put(rqc);
00889 }
00890 else
00891 {
00892 bcRqcToContext_[FunctionalOnly].put(rqc, context);
00893 bcRegionQuadComboNonzeroDerivs_[FunctionalOnly].put(rqc, nonzeros);
00894 }
00895 }
00896
00897
00898
00899
00900
00901 if (compTypes_.contains(FunctionalAndGradient))
00902 {
00903 Tabs tab3;
00904 SUNDANCE_MSG2(rqcVerb, tab3 << "preparing functional and gradient calculation");
00905 EvalContext context(rqc, makeSet(0,1), contextID[3]);
00906 context.setSetupVerbosity(symbVerb);
00907 DerivSet nonzeros;
00908 nonzeros = SymbPreprocessor
00909 ::setupGradient(term,
00910 toList(vars),
00911 toList(varLinearizationPts),
00912 fixedParams,
00913 fixedParamEvalPts,
00914 toList(fixedFields),
00915 toList(fixedFieldValues),
00916 context,
00917 FunctionalAndGradient);
00918
00919 SUNDANCE_MSG2(rqcVerb, tab3 << "nonzeros are " << nonzeros);
00920 if (nonzeros.size()==0)
00921 {
00922 bcRqcToSkip_[FunctionalAndGradient].put(rqc);
00923 }
00924 else
00925 {
00926 bcRqcToContext_[FunctionalAndGradient].put(rqc, context);
00927 bcRegionQuadComboNonzeroDerivs_[FunctionalAndGradient].put(rqc, nonzeros);
00928 }
00929 }
00930 }
00931 }
00932
00933
00934
00935
00936
00937 regionQuadCombos_ = rqcSet.elements();
00938 bcRegionQuadCombos_ = rqcBCSet.elements();
00939
00940
00941
00942
00943 }
00944
00945
00946 void EquationSet
00947 ::addToVarUnkPairs(const OrderedHandle<CellFilterStub>& domain,
00948 const Set<int>& vars,
00949 const Set<int>& unks,
00950 const DerivSet& nonzeros,
00951 bool isBC,
00952 int verb)
00953 {
00954 Tabs tab;
00955 SUNDANCE_MSG2(verb, tab << "finding var-unk pairs "
00956 "for domain " << domain);
00957 SUNDANCE_MSG2(verb, tab << "isBC=" << isBC);
00958
00959 RCP<Set<OrderedPair<int, int> > > funcPairs;
00960 Map<OrderedHandle<CellFilterStub>, RCP<Set<OrderedPair<int, int> > > >* theMap;
00961
00962 if (isBC)
00963 {
00964 theMap = &(bcVarUnkPairsOnRegions_);
00965 }
00966 else
00967 {
00968 theMap = &(varUnkPairsOnRegions_);
00969 }
00970
00971 if (theMap->containsKey(domain))
00972 {
00973 funcPairs = theMap->get(domain);
00974 }
00975 else
00976 {
00977 funcPairs = rcp(new Set<OrderedPair<int, int> >());
00978 theMap->put(domain, funcPairs);
00979 }
00980
00981 for (DerivSet::const_iterator i=nonzeros.begin(); i!=nonzeros.end(); i++)
00982 {
00983 Tabs tab1;
00984 const MultipleDeriv& md = *i;
00985 if (md.order() != 2) continue;
00986
00987 Array<Deriv> f;
00988 for (MultipleDeriv::const_iterator j=md.begin(); j != md.end(); j++)
00989 {
00990 const Deriv& d = *j;
00991 TEST_FOR_EXCEPTION(!d.isFunctionalDeriv(),
00992 InternalError, "non-functional deriv "
00993 << d << " detected in EquationSet::"
00994 "addToVarUnkPairs()");
00995 f.append(d);
00996 }
00997
00998 SUNDANCE_MSG2(verb, tab1 << "f1=" << f[0].dofID()
00999 << ", f2=" << f[1].dofID() << ", vars=" << vars
01000 << ", unks=" << unks);
01001
01002 bool gotIt=false;
01003 if (unks.contains(f[0].dofID())
01004 && vars.contains(f[1].dofID()))
01005 {
01006 int unkID = f[0].dofID();
01007 int varID = f[1].dofID();
01008 funcPairs->put(OrderedPair<int, int>(varID, unkID));
01009 gotIt=true;
01010 }
01011 if (unks.contains(f[1].dofID())
01012 && vars.contains(f[0].dofID()))
01013 {
01014 int unkID = f[1].dofID();
01015 int varID = f[0].dofID();
01016 funcPairs->put(OrderedPair<int, int>(varID, unkID));
01017 gotIt=true;
01018 }
01019 TEST_FOR_EXCEPTION(!gotIt, InternalError,
01020 "no valid (var,unk) pair could be extracted from "
01021 "derivative " << md);
01022 }
01023
01024 SUNDANCE_MSG2(verb, tab << "found " << *funcPairs);
01025
01026 }
01027
01028 bool EquationSet::hasActiveWatchFlag() const
01029 {
01030 for (int i=0; i<regionQuadCombos().size(); i++)
01031 {
01032 if (regionQuadCombos()[i].watch().isActive()) return true;
01033 }
01034 for (int i=0; i<bcRegionQuadCombos().size(); i++)
01035 {
01036 if (bcRegionQuadCombos()[i].watch().isActive()) return true;
01037 }
01038 return false;
01039 }
01040
01041
01042 int EquationSet::maxWatchFlagSetting(const std::string& param) const
01043 {
01044 int rtn = 0;
01045 if (!hasActiveWatchFlag()) return 0;
01046 for (int i=0; i<regionQuadCombos().size(); i++)
01047 {
01048 int v = regionQuadCombos()[i].watch().param(param);
01049 if (v > rtn) rtn = v;
01050 }
01051 for (int i=0; i<bcRegionQuadCombos().size(); i++)
01052 {
01053 int v = bcRegionQuadCombos()[i].watch().param(param);
01054 if (v > rtn) rtn = v;
01055 }
01056 return rtn;
01057 }
01058
01059 Array<Expr> EquationSet::flattenSpectral(const Array<Expr>& expr) const
01060 {
01061 Array<Expr> rtn(expr.size());
01062 for (int i=0; i<expr.size(); i++)
01063 {
01064 const Expr& e = expr[i];
01065 rtn[i] = flattenSpectral(e);
01066 }
01067 return rtn;
01068 }
01069
01070 Expr EquationSet::flattenSpectral(const Expr& expr) const
01071 {
01072 Array<Expr> rtn(expr.size());
01073 for (int i=0; i<expr.size(); i++)
01074 {
01075 if (expr[i].size() == 1)
01076 {
01077 const SpectralExpr* se
01078 = dynamic_cast<const SpectralExpr*>(expr[i][0].ptr().get());
01079 if (se != 0)
01080 {
01081 int nt = se->getSpectralBasis().nterms();
01082 Array<Expr> e(nt);
01083 for (int j=0; j<nt; j++)
01084 {
01085 e[j] = se->getCoeff(j);
01086 }
01087 rtn[i] = new ListExpr(e);
01088 }
01089 else
01090 {
01091 rtn[i] = expr[i];
01092 }
01093 }
01094 else
01095 {
01096 rtn[i] = flattenSpectral(expr[i]);
01097 }
01098 }
01099 Expr r = new ListExpr(rtn);
01100 return r.flatten();
01101
01102 }
01103
01104 const RCP<Set<OrderedPair<int, int> > >& EquationSet::
01105 bcVarUnkPairs(const OrderedHandle<CellFilterStub>& domain) const
01106 {
01107 TEST_FOR_EXCEPTION(!bcVarUnkPairsOnRegions_.containsKey(domain),
01108 InternalError,
01109 "equation set does not have a var-unk pair list for "
01110 "bc region " << domain);
01111 const RCP<Set<OrderedPair<int, int> > >& rtn
01112 = bcVarUnkPairsOnRegions_.get(domain);
01113
01114 TEST_FOR_EXCEPTION(rtn.get()==0, InternalError,
01115 "null var-unk pair list for BC region " << domain);
01116 return rtn;
01117 }
01118
01119 bool EquationSet::isBCRegion(int d) const
01120 {
01121 return fsr_->isBCRegion(d);
01122 }
01123
01124
01125 EvalContext EquationSet::rqcToContext(ComputationType compType,
01126 const RegionQuadCombo& r) const
01127 {
01128 TEST_FOR_EXCEPTION(!rqcToContext_.containsKey(compType),
01129 InternalError,
01130 "EquationSet::rqcToContext() did not find key "
01131 << compType);
01132 TEST_FOR_EXCEPTION(!rqcToContext_.get(compType).containsKey(r),
01133 InternalError,
01134 "EquationSet::rqcToContext(" << compType
01135 << ") did not find expected key "
01136 << r);
01137
01138 return rqcToContext_.get(compType).get(r);
01139 }
01140
01141 EvalContext EquationSet::bcRqcToContext(ComputationType compType,
01142 const RegionQuadCombo& r) const
01143 {
01144 TEST_FOR_EXCEPTION(!bcRqcToContext_.containsKey(compType),
01145 InternalError,
01146 "EquationSet::bcRqcToContext() did not find key "
01147 << compType);
01148 TEST_FOR_EXCEPTION(!bcRqcToContext_.get(compType).containsKey(r),
01149 InternalError,
01150 "EquationSet::bcRqcToContext(" << compType
01151 << ") did not find expected key "
01152 << r);
01153 return bcRqcToContext_.get(compType).get(r);
01154 }
01155
01156
01157 bool EquationSet::skipRqc(ComputationType compType,
01158 const RegionQuadCombo& r) const
01159 {
01160 TEST_FOR_EXCEPTION(!rqcToSkip_.containsKey(compType),
01161 InternalError,
01162 "EquationSet::skipRqc() did not find expected key "
01163 << compType);
01164
01165 return rqcToSkip_.get(compType).contains(r);
01166 }
01167
01168 bool EquationSet::skipBCRqc(ComputationType compType,
01169 const RegionQuadCombo& r) const
01170 {
01171 TEST_FOR_EXCEPTION(!bcRqcToSkip_.containsKey(compType),
01172 InternalError,
01173 "EquationSet::skipBCRqc() did not find expected key "
01174 << compType);
01175
01176 return bcRqcToSkip_.get(compType).contains(r);
01177 }
01178
01179 const DerivSet& EquationSet::nonzeroFunctionalDerivs(ComputationType compType,
01180 const RegionQuadCombo& r) const
01181 {
01182 TEST_FOR_EXCEPTION(!regionQuadComboNonzeroDerivs_.containsKey(compType),
01183 InternalError,
01184 "EquationSet:nonzeroFunctionalDerivs() did not find key "
01185 << compType);
01186 return regionQuadComboNonzeroDerivs_.get(compType).get(r);
01187 }
01188
01189 const DerivSet& EquationSet::nonzeroBCFunctionalDerivs(ComputationType compType,
01190 const RegionQuadCombo& r) const
01191 {
01192 TEST_FOR_EXCEPTION(!bcRegionQuadComboNonzeroDerivs_.containsKey(compType),
01193 InternalError,
01194 "EquationSet:nonzeroBCFunctionalDerivs() did not find key "
01195 << compType);
01196 return bcRegionQuadComboNonzeroDerivs_.get(compType).get(r);
01197 }
01198
01199
01200
01201 int EquationSet::reducedVarID(int varID) const
01202 {
01203 return fsr_->reducedVarID(varID);
01204 }
01205
01206 int EquationSet::reducedUnkID(int unkID) const
01207 {
01208 return fsr_->reducedUnkID(unkID);
01209 }
01210
01211
01212 int EquationSet::reducedUnkParamID(int unkParamID) const
01213 {
01214 return fsr_->reducedUnkParamID(unkParamID);
01215 }
01216
01217 int EquationSet::reducedFixedParamID(int paramID) const
01218 {
01219 return fsr_->reducedFixedParamID(paramID);
01220 }
01221
01222
01223 Expr EquationSet::toList(const Array<Expr>& e)
01224 {
01225 return new ListExpr(e);
01226 }
01227
01228
01229 int EquationSet::blockForVarID(int varID) const
01230 {
01231 return fsr_->blockForVarID(varID);
01232 }
01233
01234 int EquationSet::blockForUnkID(int unkID) const
01235 {
01236 return fsr_->blockForUnkID(unkID);
01237 }
01238
01239 const Set<OrderedHandle<CellFilterStub> >& EquationSet::regionsForTestFunc(int testID) const
01240 {
01241 return fsr_->regionsForTestFunc(testID);
01242 }
01243
01244 const Set<OrderedHandle<CellFilterStub> >& EquationSet::regionsForUnkFunc(int unkID) const
01245 {
01246 return fsr_->regionsForUnkFunc(unkID);
01247 }
01248
01249 int EquationSet::indexForRegion(const OrderedHandle<CellFilterStub>& region) const
01250 {
01251 return fsr_->indexForRegion(region);
01252 }
01253
01254 int EquationSet::numRegions() const {return fsr_->numRegions();}
01255
01256 const RCP<CellFilterStub>& EquationSet::region(int d) const
01257 {return fsr_->region(d);}
01258
01259
01260
01261
01262 int EquationSet::numVarBlocks() const
01263 {return fsr_->numVarBlocks();}
01264
01265
01266 int EquationSet::numUnkBlocks() const
01267 {return fsr_->numUnkBlocks();}
01268
01269
01270 int EquationSet::numUnkParams() const
01271 {return fsr_->numUnkParams();}
01272
01273
01274 int EquationSet::numFixedParams() const
01275 {return fsr_->numFixedParams();}
01276
01277
01278 int EquationSet::numVars(int block) const
01279 {return fsr_->numVars(block);}
01280
01281
01282 int EquationSet::numUnks(int block) const
01283 {return fsr_->numUnks(block);}
01284
01285
01286 int EquationSet::numVarIDs(int block) const
01287 {return fsr_->numVarIDs(block);}
01288
01289
01290 int EquationSet::numUnkIDs(int block) const
01291 {return fsr_->numUnkIDs(block);}
01292
01293
01294 RCP<const CommonFuncDataStub> EquationSet::varFuncData(int b, int i) const
01295 {return fsr_->varFuncData(b,i);}
01296
01297
01298
01299 RCP<const CommonFuncDataStub> EquationSet::unkFuncData(int b, int i) const
01300 {return fsr_->unkFuncData(b,i);}
01301
01302
01303 const Expr& EquationSet::unkParam(int i) const
01304 {return fsr_->unkParam(i);}
01305
01306
01307
01308 bool EquationSet::hasVarID(int fid) const
01309 {return fsr_->hasVarID(fid);}
01310
01311
01312
01313 bool EquationSet::hasUnkID(int fid) const
01314 {return fsr_->hasUnkID(fid);}
01315
01316
01317
01318 bool EquationSet::hasUnkParamID(int fid) const
01319 {return fsr_->hasUnkParamID(fid);}
01320
01321
01322
01323 bool EquationSet::hasFixedParamID(int fid) const
01324 {return fsr_->hasFixedParamID(fid);}
01325
01326
01327
01328
01329 const Set<int>& EquationSet::varsOnRegion(int d) const
01330 {return fsr_->varsOnRegion(d);}
01331
01332
01333
01334 const Set<int>& EquationSet::unksOnRegion(int d) const
01335 {return fsr_->unksOnRegion(d);}
01336
01337
01338
01339
01340 const Set<int>& EquationSet::bcVarsOnRegion(int d) const
01341 {return fsr_->bcVarsOnRegion(d);}
01342
01343
01344
01345
01346 const Set<int>& EquationSet::bcUnksOnRegion(int d) const
01347 {return fsr_->bcUnksOnRegion(d);}
01348
01349
01350
01351 const Array<Set<int> >& EquationSet::reducedVarsOnRegion(const OrderedHandle<CellFilterStub>& r) const
01352 {return fsr_->reducedVarsOnRegion(r);}
01353
01354
01355
01356 const Array<Set<int> >& EquationSet::reducedUnksOnRegion(const OrderedHandle<CellFilterStub>& r) const
01357 {return fsr_->reducedVarsOnRegion(r);}
01358
01359
01360
01361
01362 int EquationSet::unreducedVarID(int block, int reducedVarID) const
01363 {return fsr_->unreducedVarID(block, reducedVarID);}
01364
01365
01366
01367 int EquationSet::unreducedUnkID(int block, int reducedUnkID) const
01368 {return fsr_->unreducedUnkID(block, reducedUnkID);}
01369
01370
01371
01372
01373 int EquationSet::unreducedUnkParamID(int reducedUnkParamID) const
01374 {return fsr_->unreducedUnkParamID(reducedUnkParamID);}
01375
01376
01377
01378
01379 int EquationSet::unreducedFixedParamID(int reducedParamID) const
01380 {return fsr_->unreducedFixedParamID(reducedParamID);}
01381