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

Site Contact