SundanceAssembler.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 "SundanceAssembler.hpp"
00043 #include "SundanceDOFMapBuilder.hpp"
00044 #include "SundanceOut.hpp"
00045 #include "PlayaTabs.hpp"
00046 #include "SundanceCellFilter.hpp"
00047 #include "SundanceCellSet.hpp"
00048 #include "SundanceTrivialGrouper.hpp"
00049 #include "SundanceDOFMapBase.hpp"
00050 #include "SundanceEquationSet.hpp"
00051 #include "SundanceDiscreteSpace.hpp"
00052 #include "SundanceDiscreteFunction.hpp"
00053 #include "SundanceIntegralGroup.hpp"
00054 #include "SundanceGrouperBase.hpp"
00055 #include "SundanceEvalManager.hpp"
00056 #include "SundanceStdFwkEvalMediator.hpp"
00057 #include "SundanceEvaluatableExpr.hpp"
00058 #include "PlayaLoadableVector.hpp"
00059 #include "PlayaLoadableMatrix.hpp"
00060 #include "SundanceQuadratureEvalMediator.hpp"
00061 #include "SundanceCurveEvalMediator.hpp"
00062 #include "SundanceEvaluator.hpp"
00063 #include "Teuchos_Time.hpp"
00064 #include "Teuchos_TimeMonitor.hpp"
00065 #include "Epetra_HashTable.h"
00066 #include "SundanceIntHashSet.hpp"
00067 #include "PlayaDefaultBlockVectorSpaceDecl.hpp"
00068 #include "PlayaBlockOperatorBaseDecl.hpp"
00069 #include "PlayaSimpleBlockOpDecl.hpp"
00070 #include "SundanceAssemblyKernelBase.hpp"
00071 #include "SundanceVectorAssemblyKernel.hpp"
00072 #include "SundanceMatrixVectorAssemblyKernel.hpp"
00073 #include "SundanceFunctionalAssemblyKernel.hpp"
00074 #include "SundanceFunctionalGradientAssemblyKernel.hpp"
00075 #include "SundanceAssemblyTransformationBuilder.hpp"
00076 #ifndef HAVE_TEUCHOS_EXPLICIT_INSTANTIATION
00077 #include "PlayaLinearOperatorImpl.hpp"
00078 #include "PlayaSimpleBlockOpImpl.hpp"
00079 #include "PlayaDefaultBlockVectorImpl.hpp"
00080 #endif
00081 
00082 
00083 
00084 using namespace Teuchos;
00085 using namespace Playa;
00086 using std::max;
00087 using std::min;
00088 
00089 
00090 static Time& assemblerCtorTimer() 
00091 {
00092   static RCP<Time> rtn 
00093     = TimeMonitor::getNewTimer("assembler ctor"); 
00094   return *rtn;
00095 }
00096 
00097 
00098 
00099 
00100 static Time& graphBuildTimer() 
00101 {
00102   static RCP<Time> rtn 
00103     = TimeMonitor::getNewTimer("matrix graph determination"); 
00104   return *rtn;
00105 }
00106 
00107 static Time& colSearchTimer() 
00108 {
00109   static RCP<Time> rtn 
00110     = TimeMonitor::getNewTimer("graph column processing"); 
00111   return *rtn;
00112 }
00113 
00114 
00115 
00116 static Time& matAllocTimer() 
00117 {
00118   static RCP<Time> rtn 
00119     = TimeMonitor::getNewTimer("matrix allocation"); 
00120   return *rtn;
00121 }
00122 
00123 static Time& matFinalizeTimer() 
00124 {
00125   static RCP<Time> rtn 
00126     = TimeMonitor::getNewTimer("matrix graph packing"); 
00127   return *rtn;
00128 }
00129 
00130 static Time& graphFlatteningTimer() 
00131 {
00132   static RCP<Time> rtn 
00133     = TimeMonitor::getNewTimer("tmp graph flattening"); 
00134   return *rtn;
00135 }
00136 
00137 
00138 
00139 Assembler
00140 ::Assembler(const Mesh& mesh, 
00141   const RCP<EquationSet>& eqn,
00142   const Array<VectorType<double> >& rowVectorType,
00143   const Array<VectorType<double> >& colVectorType,
00144   bool partitionBCs)
00145   : partitionBCs_(partitionBCs),
00146     numConfiguredColumns_(0),
00147     mesh_(mesh),
00148     eqn_(eqn),
00149     rowMap_(),
00150     colMap_(),
00151     externalRowSpace_(eqn->numVarBlocks()),
00152     externalColSpace_(eqn->numUnkBlocks()),
00153     privateRowSpace_(eqn->numVarBlocks()),
00154     privateColSpace_(eqn->numUnkBlocks()),
00155     bcRows_(eqn->numVarBlocks()),
00156     rqc_(),
00157     contexts_(),
00158     isBCRqc_(),
00159     isInternalBdry_(),
00160     groups_(),
00161     mediators_(),
00162     evalExprs_(),
00163     evalMgr_(rcp(new EvalManager())),
00164     isBCRow_(eqn->numVarBlocks()),
00165     isBCCol_(eqn->numUnkBlocks()),
00166     remoteBCCols_(eqn->numUnkBlocks()),
00167     lowestRow_(eqn->numVarBlocks()),
00168     lowestCol_(eqn->numUnkBlocks()),
00169     rowVecType_(rowVectorType),
00170     colVecType_(colVectorType),
00171     testIDToBlockMap_(),
00172     unkIDToBlockMap_()
00173 {
00174   TimeMonitor timer(assemblerCtorTimer());
00175   init(mesh, eqn);
00176 }
00177 
00178 Assembler
00179 ::Assembler(const Mesh& mesh, 
00180   const RCP<EquationSet>& eqn)
00181   : partitionBCs_(false),
00182     numConfiguredColumns_(0),
00183     mesh_(mesh),
00184     eqn_(eqn),
00185     rowMap_(),
00186     colMap_(),
00187     externalRowSpace_(eqn->numVarBlocks()),
00188     externalColSpace_(eqn->numUnkBlocks()),
00189     privateRowSpace_(eqn->numVarBlocks()),
00190     privateColSpace_(eqn->numUnkBlocks()),
00191     bcRows_(eqn->numVarBlocks()),
00192     rqc_(),
00193     contexts_(),
00194     isBCRqc_(),
00195     isInternalBdry_(),
00196     groups_(),
00197     mediators_(),
00198     evalExprs_(),
00199     evalMgr_(rcp(new EvalManager())),
00200     isBCRow_(eqn->numVarBlocks()),
00201     isBCCol_(eqn->numUnkBlocks()),
00202     remoteBCCols_(eqn->numUnkBlocks()),
00203     lowestRow_(eqn->numVarBlocks()),
00204     lowestCol_(eqn->numUnkBlocks()),
00205     rowVecType_(),
00206     colVecType_(),
00207     testIDToBlockMap_(),
00208     unkIDToBlockMap_(),
00209     fixedParamIDToVectorNumber_()
00210 {
00211   TimeMonitor timer(assemblerCtorTimer());
00212   init(mesh, eqn);
00213 }
00214 
00215 // Utility function
00216 namespace {
00217 const Set<int> &
00218 findNonZeroIndices(const Array<int> &source, Set<int> &result)
00219 {
00220   result.clear();
00221   for (int i=0; i<source.size(); ++i)
00222   {
00223     if (source[i] != 0)
00224     {
00225       result.insert(result.end(), i);
00226     }
00227   }
00228   return result;
00229 }
00230 }
00231 
00232 void Assembler::init(const Mesh& mesh, const RCP<EquationSet>& eqn)
00233 {
00234   Tabs tab0(0);
00235 
00236   /* Decide a verbosity level for the overall setup */
00237   int verb = 0;
00238   if (eqn->hasActiveWatchFlag())
00239     verb = eqn->maxWatchFlagSetting("assembler setup");
00240 
00241   SUNDANCE_BANNER1(verb, tab0, "Assembler setup");
00242 
00243 
00244   /* Create an integral grouper */
00245   RCP<GrouperBase> grouper 
00246     = rcp(new TrivialGrouper());
00247 
00248 
00249   /* Find out which types of computations this assembler will 
00250    * be required to carry out */
00251   const Set<ComputationType>& compTypes = eqn->computationTypes();
00252 
00253 
00254   /* Create the DOF map for the row space */
00255   DOFMapBuilder mapBuilder(eqn->maxWatchFlagSetting("dof map setup"));
00256 
00257   if (compTypes.contains(VectorOnly) 
00258     || compTypes.contains(Sensitivities) 
00259     || compTypes.contains(FunctionalAndGradient))
00260   {
00261     Tabs tab1;
00262     SUNDANCE_MSG2(verb, tab1 << "building row spaces");
00263     mapBuilder = DOFMapBuilder(mesh, eqn->fsr(), partitionBCs_, 
00264       eqn->maxWatchFlagSetting("dof map setup"));
00265 
00266     rowMap_ = mapBuilder.rowMap();
00267     isBCRow_ = mapBuilder.isBCRow();
00268     isBCCol_ = mapBuilder.isBCCol();
00269 
00270     bcRows_.resize(isBCRow_.size());
00271     for (int b=0; b<isBCRow_.size(); ++b)
00272     {
00273       bcRows_[b] = rcp(new Set<int>);
00274       findNonZeroIndices(*isBCRow_[b], *bcRows_[b]);
00275     }
00276 
00277     lowestRow_.resize(eqn_->numVarBlocks());
00278     /* create discrete space for each block */
00279     for (int b=0; b<eqn_->numVarBlocks(); b++) 
00280     {
00281       Tabs tab2;
00282       lowestRow_[b] = rowMap_[b]->lowestLocalDOF();
00283       SUNDANCE_MSG2(verb, tab2 << "block " << b << ": lowest row="
00284         << lowestRow_[b]);
00285       externalRowSpace_[b] = rcp(
00286         new DiscreteSpace(mesh, testBasisArray(mapBuilder.fsr())[b], 
00287           rowMap_[b], rowVecType_[b]));
00288       privateRowSpace_[b] = externalRowSpace_[b];
00289       SUNDANCE_MSG2(verb, tab2 << "block " << b << ": done forming row space");
00290     }
00291   }
00292 
00293   if (!eqn->isFunctionalCalculator())
00294   {
00295     Tabs tab1;
00296     /* Create the DOF map for the column space */
00297     SUNDANCE_MSG2(verb, tab1 << "building column spaces");
00298     colMap_ = mapBuilder.colMap();
00299     /* create discrete space for each block */
00300     for (int b=0; b<eqn_->numUnkBlocks(); b++) 
00301     {
00302       Tabs tab2;
00303       externalColSpace_[b] 
00304         = rcp(new DiscreteSpace(mesh, unkBasisArray(mapBuilder.fsr())[b], 
00305             colMap_[b], colVecType_[b]));
00306       privateColSpace_[b] = externalColSpace_[b];
00307       SUNDANCE_MSG2(verb, tab2 << "block " << b << ": done forming col space");
00308     }
00309 
00310     /* initialize empty tables of information for each RQC in a 
00311      * matrix-vector calculation */
00312     groups_.put(MatrixAndVector, Array<Array<RCP<IntegralGroup> > >());
00313     rqcRequiresMaximalCofacets_.put(MatrixAndVector, 
00314       Array<IntegrationCellSpecifier>());
00315     skipRqc_.put(MatrixAndVector, Array<int>());
00316     contexts_.put(MatrixAndVector, Array<EvalContext>());
00317     evalExprs_.put(MatrixAndVector, Array<const EvaluatableExpr*>());
00318 
00319     /* create tables for vector calculation */
00320     groups_.put(VectorOnly, Array<Array<RCP<IntegralGroup> > >());
00321     rqcRequiresMaximalCofacets_.put(VectorOnly, 
00322       Array<IntegrationCellSpecifier>());
00323     skipRqc_.put(VectorOnly, Array<int>());
00324     contexts_.put(VectorOnly, Array<EvalContext>());
00325     evalExprs_.put(VectorOnly, Array<const EvaluatableExpr*>());
00326 
00327     if (eqn->isSensitivityCalculator())
00328     {
00329       fixedParamIDToVectorNumber_ 
00330         = eqn->fsr()->fixedParamIDToReducedFixedParamIDMap();
00331 
00332       /* create tables for sensitivity calculation */
00333       groups_.put(Sensitivities, Array<Array<RCP<IntegralGroup> > >());
00334       rqcRequiresMaximalCofacets_.put(Sensitivities, 
00335         Array<IntegrationCellSpecifier>());
00336       skipRqc_.put(Sensitivities, Array<int>());
00337       contexts_.put(Sensitivities, Array<EvalContext>());
00338       evalExprs_.put(Sensitivities, Array<const EvaluatableExpr*>());
00339       
00340     }
00341   }
00342   else
00343   {
00344     /* create tables for functional and gradient calculation */
00345     groups_.put(FunctionalAndGradient, Array<Array<RCP<IntegralGroup> > >());
00346     rqcRequiresMaximalCofacets_.put(FunctionalAndGradient, 
00347       Array<IntegrationCellSpecifier>());
00348     skipRqc_.put(FunctionalAndGradient, Array<int>());
00349     contexts_.put(FunctionalAndGradient, Array<EvalContext>());
00350     evalExprs_.put(FunctionalAndGradient, Array<const EvaluatableExpr*>());
00351     /* create tables for functional calculation */
00352     groups_.put(FunctionalOnly, Array<Array<RCP<IntegralGroup> > >());
00353     rqcRequiresMaximalCofacets_.put(FunctionalOnly, 
00354       Array<IntegrationCellSpecifier>());
00355     skipRqc_.put(FunctionalOnly, Array<int>());
00356     contexts_.put(FunctionalOnly, Array<EvalContext>());
00357     evalExprs_.put(FunctionalOnly, Array<const EvaluatableExpr*>());
00358   }
00359 
00360 
00361 
00362 
00363 
00364   /* --- We now loop over non-BC RQCs, doing initialization tasks for each */
00365   SUNDANCE_MSG1(verb, tab0 << std::endl 
00366     << tab0 << "=== setting up non-BC region-quadrature combinations");
00367 
00368   for (int r=0; r<eqn->regionQuadCombos().size(); r++)
00369   {
00370     Tabs tab1;
00371     Tabs tab12;
00372     const RegionQuadCombo& rqc = eqn->regionQuadCombos()[r];
00373 
00374     /* Determine the verbosity setting for this RQC */
00375     bool watchMe = rqc.watch().isActive();
00376 
00377     int rqcVerb = verb;
00378     int integralCtorVerb = 0;
00379     int integrationVerb = 0;
00380     int integralTransformVerb = 0;
00381     if (watchMe) 
00382     {
00383       rqcVerb = max(4,rqcVerb);
00384       integralCtorVerb = rqc.watch().param("integration setup");
00385       integrationVerb = rqc.watch().param("integration");
00386       integralTransformVerb = rqc.watch().param("integral transformation");
00387     }
00388 
00389 
00390     /* Note that I'm not an essential BC */
00391     rqc_.append(rqc);
00392     isBCRqc_.append(false);
00393 
00394     /* Find the expression for this RQC */
00395     const Expr& expr = eqn->expr(rqc);
00396 
00397     SUNDANCE_MSG2(rqcVerb, tab1 << std::endl << tab1 << "------------------------------------------------");
00398     SUNDANCE_MSG2(rqcVerb, tab1 << "initializing assembly for"
00399       << std::endl << tab12 << "r=" << r << " rqc=" 
00400       << rqc << std::endl << tab12 << std::endl << tab12 << "------------------------------"
00401       << std::endl << tab12 << "expr = " << expr
00402       << std::endl << tab12 << "------------------------------"
00403       );
00404 
00405     
00406     /* Find the cell type needed for this RQC */
00407     int cellDim = CellFilter(rqc.domain()).dimension(mesh);
00408     CellType cellType = mesh.cellType(cellDim);
00409     CellType maxCellType = mesh.cellType(mesh.spatialDim());
00410     QuadratureFamily quad(rqc.quad());
00411 
00412     /* Detect internal boundaries. These need special handling */
00413     bool isInternalBdry = detectInternalBdry(cellDim, rqc.domain());
00414     isInternalBdry_.append(isInternalBdry);
00415 
00416     SUNDANCE_MSG2(rqcVerb, tab12 << "isInternalBdry=" << isInternalBdry);
00417 
00418     /* Do setup for each required computation type */
00419     bool rqcUsed = false;
00420 
00421     for (Set<ComputationType>::const_iterator 
00422            i=eqn->computationTypes().begin(); 
00423          i!=eqn->computationTypes().end();
00424          i++)
00425     {
00426       Tabs tab2;
00427       const ComputationType& compType = *i;
00428       SUNDANCE_MSG2(rqcVerb, tab12 << std::endl << tab12
00429         << "** computation type " << compType);
00430       /* Some RQCs may be unused in a given calculation. For example, an RQC
00431        * may be needed for vector calculation but not matrix-vector 
00432        * calculation. See if this RQC is needed for the present 
00433        * computation type. If not, there's nothing more to do here. */
00434       if (eqn->skipRqc(compType, rqc))
00435       {
00436         SUNDANCE_MSG2(rqcVerb, tab2 << "RQC not needed for computation type  " 
00437           << compType);
00438         skipRqc_[compType].append(true);
00439         EvalContext dummy;
00440         const EvaluatableExpr* dummyEx = 0;
00441         Array<RCP<IntegralGroup> > dummyGroups;
00442         IntegrationCellSpecifier dummyCellSpec ;
00443         contexts_[compType].append(dummy);
00444         evalExprs_[compType].append(dummyEx);
00445         groups_[compType].append(dummyGroups);
00446         rqcRequiresMaximalCofacets_[compType].append(dummyCellSpec);
00447       }
00448       else
00449       {
00450         /* If we're to this point, we know the RQC is needed for this 
00451          * computation type */
00452         rqcUsed = true;
00453         skipRqc_[compType].append(false);
00454 
00455         /* Look up a "context" object that we'll use as a key for different
00456          * evaluations of this expression */
00457         EvalContext context = eqn->rqcToContext(compType, rqc);
00458 
00459         SUNDANCE_MSG2(rqcVerb, tab2 << "context " << context.brief());
00460 
00461         /* Register the context */
00462         contexts_[compType].append(context);
00463 
00464         /* Register the expression */
00465         const EvaluatableExpr* ee = EvaluatableExpr::getEvalExpr(expr);
00466         evalExprs_[compType].append(ee);
00467 
00468         /* Get the "sparsity superset" which is a description of all 
00469          * derivatives that must be computed by this expression in the
00470          * present context */
00471         const RCP<SparsitySuperset>& sparsity 
00472           = ee->sparsitySuperset(context);
00473         SUNDANCE_MSG3(rqcVerb, tab2 << "sparsity pattern " << *sparsity);
00474 
00475         /* We're now ready to create integration groups for doing the 
00476          * integrals needed in this computation for the present RQC. */
00477         Array<RCP<IntegralGroup> > groups;
00478         grouper->setVerb(integralCtorVerb, integrationVerb, integralTransformVerb);
00479         grouper->findGroups(*eqn, maxCellType, mesh.spatialDim(),
00480           cellType, cellDim, quad, sparsity, isInternalBdry, groups , rqc.paramCurve() , mesh_ );
00481         grouper->setVerb(0,0,0);
00482         groups_[compType].append(groups);
00483 
00484         /* Record whether or not integrations need to be done by reference
00485          * to maximal cofacets. */ 
00486         IntegrationCellSpecifier cellSpec 
00487           = whetherToUseCofacets(groups, ee, 
00488             cellDim==mesh_.spatialDim(), rqcVerb);
00489         SUNDANCE_MSG2(rqcVerb, tab2 << "integration: " << cellSpec);
00490         rqcRequiresMaximalCofacets_[compType].append(cellSpec);
00491       }
00492       SUNDANCE_MSG2(rqcVerb, tab12
00493         << "done with computation type " << compType);
00494     }
00495     
00496     /* If this RQC has never been used, we've made a mistake */
00497     /* Actually, no! Some terms might be unused in reduced-space methods */
00498 //    TEUCHOS_TEST_FOR_EXCEPTION(!rqcUsed, std::logic_error, "rqc=" << rqc 
00499 //      << " never used for any computation???");
00500 
00501     if (rqcUsed)
00502     {
00503       SUNDANCE_MSG2(rqcVerb, tab12 << "creating evaluation mediator for rqc=" 
00504         << rqc);
00505       SUNDANCE_MSG2(rqcVerb, tab12 << "expr = " << expr);
00506       
00507       /* Finally, create an evaluation mediator for this RQC. The evaluation
00508        * mediator is the object through which symbolic objects refer to
00509        * mesh-dependent quantities (e.g., discrete functions) during
00510        * evaluation.  , if we have to evaluate a curve integral then
00511        * use a special mediator */
00512       if ( rqc.paramCurve().isCurveIntegral() && rqc.paramCurve().isCurveValid() )
00513       { // ----- curve Integral ------
00514          mediators_.append(rcp(new CurveEvalMediator(mesh, rqc.paramCurve() , cellDim, quad)));
00515       }
00516       else
00517       {
00518          mediators_.append(rcp(new QuadratureEvalMediator(mesh, cellDim, quad)));
00519       }
00520     }
00521     else
00522     {
00523         SUNDANCE_MSG2(rqcVerb, tab1 << "creating empty eval mediator for unused rqc");
00524         mediators_.append(RCP<StdFwkEvalMediator>());
00525     }
00526     SUNDANCE_MSG2(rqcVerb, tab1 
00527       << "done with RQC");
00528   }
00529 
00530 
00531 
00532   /* --- We now loop over BC RQCs, doing initialization tasks for each */
00533   SUNDANCE_MSG1(verb, tab0 << std::endl 
00534     << tab0 << "=== setting up BC region-quadrature combinations");
00535   
00536   for (int r=0; r<eqn->bcRegionQuadCombos().size(); r++)
00537   {
00538     Tabs tab1;
00539     const RegionQuadCombo& rqc = eqn->bcRegionQuadCombos()[r];
00540 
00541     /* Determine the verbosity setting for this RQC */
00542     bool watchMe = rqc.watch().isActive();
00543     int rqcVerb = verb;
00544     int integralCtorVerb = 0;
00545     int integrationVerb = 0;
00546     int integralTransformVerb = 0;
00547     if (watchMe) 
00548     {
00549       rqcVerb = max(4,rqcVerb);
00550       integralCtorVerb = rqc.watch().param("integration setup");
00551       integrationVerb = rqc.watch().param("integration");
00552       integralTransformVerb = rqc.watch().param("integral transformation");
00553     }
00554 
00555     /* Note that I am an essential BC */
00556     rqc_.append(rqc);
00557     isBCRqc_.append(true);
00558 
00559 
00560     /* Find the expression for this RQC */    
00561     const Expr& expr = eqn->bcExpr(rqc);
00562 
00563     SUNDANCE_MSG2(rqcVerb, tab1 << std::endl << tab1 
00564       << "------------------------------------------------");
00565     SUNDANCE_MSG1(rqcVerb, tab1 << "initializing assembly for BC r=" << r
00566       << " rqc=" 
00567       << rqc << std::endl << tab1 
00568       << "expr = " << expr);
00569       
00570     /* Find the cell type needed for this RQC */    
00571     int cellDim = CellFilter(rqc.domain()).dimension(mesh);
00572     CellType cellType = mesh.cellType(cellDim);
00573     CellType maxCellType = mesh.cellType(mesh.spatialDim());
00574     QuadratureFamily quad(rqc.quad());
00575 
00576     /* Detect internal boundaries. These need special handling */
00577     bool isInternalBdry = detectInternalBdry(cellDim, rqc.domain());
00578     isInternalBdry_.append(isInternalBdry);
00579 
00580     /* Do setup for each required computation type */
00581     bool rqcUsed = false;
00582 
00583     for (Set<ComputationType>::const_iterator 
00584            i=eqn->computationTypes().begin(); 
00585          i!=eqn->computationTypes().end();
00586          i++)
00587     {
00588       Tabs tab2;
00589       const ComputationType& compType = *i;
00590       SUNDANCE_MSG2(rqcVerb, tab1 << std::endl << tab1 
00591         << "** computation type " << compType);
00592       
00593       /* Some RQCs may be unused in a given calculation. For example, an RQC
00594        * may be needed for vector calculation but not matrix-vector 
00595        * calculation. See if this RQC is needed for the present 
00596        * computation type. If not, there's nothing more to do here. */
00597       if (eqn->skipBCRqc(compType, rqc))
00598       {
00599         SUNDANCE_MSG2(rqcVerb, 
00600           tab2 << "this rqc not needed for computation type " << compType);
00601         skipRqc_[compType].append(true);
00602         EvalContext dummy;
00603         const EvaluatableExpr* dummyEx = 0;
00604         Array<RCP<IntegralGroup> > dummyGroups;
00605         IntegrationCellSpecifier dummyCellSpec ;
00606         contexts_[compType].append(dummy);
00607         evalExprs_[compType].append(dummyEx);
00608         groups_[compType].append(dummyGroups);
00609         rqcRequiresMaximalCofacets_[compType].append(dummyCellSpec);
00610       }
00611       else
00612       {
00613         /* If we're to this point, we know the RQC is needed for this 
00614          * computation type */   
00615         rqcUsed = true;
00616         skipRqc_[compType].append(false);
00617 
00618         /* Look up a "context" object that we'll use as a key for different
00619          * evaluations of this expression */
00620         EvalContext context = eqn->bcRqcToContext(compType, rqc);
00621         SUNDANCE_MSG2(rqcVerb, tab2 << "context " << context);
00622 
00623       
00624         contexts_[compType].append(context);
00625         const EvaluatableExpr* ee = EvaluatableExpr::getEvalExpr(expr);
00626         evalExprs_[compType].append(ee);
00627         const RCP<SparsitySuperset>& sparsity 
00628           = ee->sparsitySuperset(context);
00629         SUNDANCE_MSG3(rqcVerb, tab2 << "sparsity pattern " << *sparsity);
00630 
00631         Array<RCP<IntegralGroup> > groups;
00632         grouper->setVerb(integralCtorVerb, integrationVerb, integralTransformVerb);
00633         grouper->findGroups(*eqn, maxCellType, mesh.spatialDim(),
00634           cellType, cellDim, quad, sparsity, isInternalBdry, groups , rqc.paramCurve() , mesh_ );
00635         grouper->setVerb(0,0,0);
00636         groups_[compType].append(groups);
00637         IntegrationCellSpecifier cellSpec 
00638           = whetherToUseCofacets(groups, ee, 
00639             cellDim==mesh_.spatialDim(), rqcVerb);
00640         SUNDANCE_MSG2(rqcVerb, tab2 << "integration: " << cellSpec);
00641         rqcRequiresMaximalCofacets_[compType].append(cellSpec);
00642         SUNDANCE_MSG2(rqcVerb, tab1
00643           << "done with computation type " << compType);
00644       }
00645 /* Turn this test off, because some terms might be unused in reduced-space
00646  * methods */
00647 //      TEUCHOS_TEST_FOR_EXCEPTION(!rqcUsed, std::logic_error, "BC rqc=" << rqc 
00648 //        << " never used for any computation???");
00649       if (rqcUsed)
00650       {
00651         SUNDANCE_MSG2(rqcVerb, tab1 << "creating evaluation mediator for BC rqc=" 
00652           << rqc << std::endl << tab1 << "expr = " << expr);
00653         // in case of curve integral use a special mediator
00654         if ( rqc.paramCurve().isCurveIntegral() && rqc.paramCurve().isCurveValid() )
00655         { // ----- curve Integral ------
00656           mediators_.append(rcp(new CurveEvalMediator(mesh, rqc.paramCurve(), cellDim, quad)));
00657         }
00658         else
00659         {
00660           mediators_.append(rcp(new QuadratureEvalMediator(mesh, cellDim, quad)));
00661         }
00662       }
00663       else
00664       {
00665         SUNDANCE_MSG2(rqcVerb, tab1 << "creating empty eval mediator for unused BC rqc");
00666         mediators_.append(RCP<StdFwkEvalMediator>());
00667       }
00668     }
00669     SUNDANCE_MSG2(rqcVerb, tab1 
00670       << "done with BC RQC");
00671   }
00672 }
00673 
00674 bool Assembler::detectInternalBdry(int cellDim,
00675   const CellFilter& filter) const
00676 {
00677   int d = mesh_.spatialDim();
00678   if (cellDim == d-1)
00679   {
00680     CellSet cells = filter.getCells(mesh_);
00681     for (CellIterator c=cells.begin(); c!=cells.end(); c++)
00682     {
00683       if (mesh_.numMaxCofacets(cellDim, *c) > 1) return true;
00684     }      
00685   }
00686   return false;
00687 }
00688 
00689 IntegrationCellSpecifier Assembler::whetherToUseCofacets(
00690   const Array<RCP<IntegralGroup> >& groups,
00691   const EvaluatableExpr* ee,
00692   bool isMaximalCell,
00693   int verb) const
00694 {
00695   Tabs tab;
00696   SUNDANCE_MSG2(verb, 
00697     tab << "deciding whether to use cofacet cells for some integrations");
00698 
00699   if (isMaximalCell)
00700   {
00701     Tabs tab1;
00702     SUNDANCE_MSG2(verb, 
00703       tab1 << "cofacets not needed because cells are maximal");
00704     return NoTermsNeedCofacets;
00705   }
00706   
00707   IntegrationCellSpecifier cellSpec = SomeTermsNeedCofacets;
00708 
00709   bool allTermsNeedCofacets = true;
00710   bool noTermsNeedCofacets = true;
00711   for (int g=0; g<groups.size(); g++)
00712   {
00713     Tabs tab1;
00714     switch(groups[g]->usesMaximalCofacets()) 
00715     {
00716       case NoTermsNeedCofacets:
00717         allTermsNeedCofacets = false;
00718         SUNDANCE_MSG2(verb, 
00719           tab1 << "integral group " << g << " does not need cofacets");
00720         break;
00721       case AllTermsNeedCofacets:
00722       case SomeTermsNeedCofacets:
00723         noTermsNeedCofacets = false;
00724         SUNDANCE_MSG2(verb, 
00725           tab1 << "integral group " << g << " needs cofacets");
00726         break;
00727       default:
00728         TEUCHOS_TEST_FOR_EXCEPT(1);
00729     }
00730   } 
00731 
00732   if (allTermsNeedCofacets)
00733   {
00734     cellSpec = AllTermsNeedCofacets;
00735   }
00736   else if (noTermsNeedCofacets)
00737   {
00738     cellSpec = NoTermsNeedCofacets;
00739   }
00740 
00741   if (!isMaximalCell && ee->maxDiffOrderOnDiscreteFunctions() > 0)
00742   {
00743     Tabs tab1;
00744     SUNDANCE_MSG2(verb, tab1 
00745       << "(*) discrete functions will require cofacet-based transformations");
00746     if (cellSpec==NoTermsNeedCofacets) 
00747     {
00748       cellSpec = SomeTermsNeedCofacets;
00749     }
00750   }
00751   
00752   SUNDANCE_MSG2(verb, tab << "found: " << cellSpec);
00753   
00754   return cellSpec;
00755 }
00756   
00757 
00758 void Assembler::configureVector(Array<Vector<double> >& b) const 
00759 {
00760   /* Start timer, stopped upon dtor */
00761   TimeMonitor timer(configTimer());
00762 
00763   Tabs tab0;
00764   int verb = eqn_->maxWatchFlagSetting("vector config");
00765   
00766   SUNDANCE_MSG1(verb, tab0 << "in Assembler::configureVector()");
00767 
00768   /* Get the vector spaces for each block of equations */
00769   Array<VectorSpace<double> > vs(eqn_->numVarBlocks());
00770   for (int i=0; i<eqn_->numVarBlocks(); i++)
00771   {
00772     vs[i] = privateRowSpace_[i]->vecSpace();
00773   }
00774   VectorSpace<double> rowSpace;
00775   
00776   /* If we have more than one block, we make a Cartesian product space containing
00777    * the spaces for each block */
00778   if (eqn_->numVarBlocks() > 1)
00779   {
00780     rowSpace = Playa::blockSpace(vs);
00781   }
00782   else /* Otherwise we have a single, monolithic vector space */
00783   {
00784     rowSpace = vs[0];
00785   }
00786 
00787   /* Create each vector in the multivector */
00788   for (int i=0; i<b.size(); i++)
00789   {
00790     /* Create the vector. Recall that the vector space is a factory used to 
00791      * create a vector of specified size and distribution */
00792     b[i] = rowSpace.createMember();
00793 
00794     /* If the vector is blocked, configure the blocks */
00795     if (!partitionBCs_ && eqn_->numVarBlocks() > 1)
00796     {
00797       /* configure the blocks */
00798       Vector<double> vecBlock;
00799       for (int br=0; br<eqn_->numVarBlocks(); br++)
00800       {
00801         configureVectorBlock(br, vecBlock);
00802         b[i].setBlock(br, vecBlock);
00803       }
00804     }
00805     else  
00806     {
00807       /* nothing to do here except check that the vector is loadable */
00808       if (!partitionBCs_)
00809       {
00810         Playa::LoadableVector<double>* lv 
00811           = dynamic_cast<Playa::LoadableVector<double>* >(b[i].ptr().get());
00812         
00813         TEUCHOS_TEST_FOR_EXCEPTION(lv == 0, std::runtime_error,
00814           "vector is not loadable in Assembler::configureVector()");
00815       }
00816       else
00817       {
00818       }
00819     }
00820   }
00821   numConfiguredColumns_ = b.size();
00822 }
00823 
00824 void Assembler::configureVectorBlock(int br, Vector<double>& b) const 
00825 {
00826   Tabs tab0;
00827   int verb = eqn_->maxWatchFlagSetting("vector config");
00828   SUNDANCE_MSG2(verb, tab0 << "in Assembler::configureVectorBlock()");
00829   VectorSpace<double> vecSpace = privateRowSpace_[br]->vecSpace();
00830 
00831   b = vecSpace.createMember();
00832   
00833   if (!partitionBCs_)
00834   {
00835     Playa::LoadableVector<double>* lv 
00836       = dynamic_cast<Playa::LoadableVector<double>* >(b.ptr().get());
00837     
00838     TEUCHOS_TEST_FOR_EXCEPTION(lv == 0, std::runtime_error,
00839       "vector block is not loadable "
00840       "in Assembler::configureVectorBlock()");
00841   }
00842 }
00843 
00844 
00845 void Assembler::configureMatrix(LinearOperator<double>& A,
00846   Array<Vector<double> >& b) const
00847 {
00848   TimeMonitor timer(configTimer());
00849   int verb = eqn_->maxWatchFlagSetting("matrix config");
00850   Tabs tab;
00851   SUNDANCE_MSG1(verb, tab << "in Assembler::configureMatrix()");
00852   
00853 
00854   SUNDANCE_MSG1(verb,  tab << "before config: A=" << A.description());
00855   if (matNeedsConfiguration())
00856   {
00857     Tabs tab0;
00858 
00859     int nRowBlocks = rowMap_.size();
00860     int nColBlocks = colMap_.size();
00861     Array<Array<int> > isNonzero = findNonzeroBlocks();
00862 
00863     if (nRowBlocks==1 && nColBlocks==1)
00864     {
00865       configureMatrixBlock(0,0,A);
00866     }
00867     else
00868     {
00869       A = makeBlockOperator(solnVecSpace(), rowVecSpace());
00870       for (int br=0; br<nRowBlocks; br++)
00871       {
00872         for (int bc=0; bc<nColBlocks; bc++)
00873         {
00874           if (isNonzero[br][bc])
00875           {
00876             LinearOperator<double> matBlock;
00877             configureMatrixBlock(br, bc, matBlock);
00878             A.setBlock(br, bc, matBlock);
00879           }
00880         }
00881       }
00882       A.endBlockFill();
00883     }
00884     cachedAssembledMatrix_ = A;
00885   }
00886   else
00887   {
00888     Tabs tab0;
00889     SUNDANCE_MSG1(verb,
00890       tab0 << "Assembler::configureMatrix() not needed, proceeding to configure vector");
00891     A = cachedAssembledMatrix_;
00892   }
00893   SUNDANCE_MSG1(verb,  tab << "after config: A=" << A.description());
00894   configureVector(b);
00895 }
00896 
00897 void Assembler::configureMatrixBlock(int br, int bc,
00898   LinearOperator<double>& A) const 
00899 {
00900   Tabs tab;
00901   TimeMonitor timer(configTimer());
00902   int verb = eqn_->maxWatchFlagSetting("matrix config");
00903 
00904   SUNDANCE_MSG1(verb, tab << "in configureMatrixBlock()");
00905   
00906   SUNDANCE_MSG2(verb, tab << "num rows = " << rowMap()[br]->numDOFs());
00907   
00908   SUNDANCE_MSG2(verb, tab << "num cols = " << colMap()[bc]->numDOFs());
00909   
00910   VectorSpace<double> rowSpace = privateRowSpace_[br]->vecSpace();
00911   VectorSpace<double> colSpace = privateColSpace_[bc]->vecSpace();
00912 
00913   RCP<MatrixFactory<double> > matFactory ;
00914 
00915   matFactory = rowVecType_[br].createMatrixFactory(colSpace, rowSpace);
00916 
00917   IncrementallyConfigurableMatrixFactory* icmf 
00918     = dynamic_cast<IncrementallyConfigurableMatrixFactory*>(matFactory.get());
00919 
00920   CollectivelyConfigurableMatrixFactory* ccmf 
00921     = dynamic_cast<CollectivelyConfigurableMatrixFactory*>(matFactory.get());
00922 
00923   TEUCHOS_TEST_FOR_EXCEPTION(ccmf==0 && icmf==0, std::runtime_error,
00924     "Neither incremental nor collective matrix structuring "
00925     "appears to be available");
00926 
00927 
00928   /* If collective structuring is the user preference, or if incremental
00929    * structuring is not supported, do collective structuring */
00930   if (false /* (icmf==0 || !matrixEliminatesRepeatedCols()) && ccmf != 0 */)
00931   {
00932     Tabs tab1;
00933     SUNDANCE_MSG2(verb, tab1 << "Assembler: doing collective matrix structuring...");
00934     Array<int> graphData;
00935     Array<int> nnzPerRow;
00936     Array<int> rowPtrs;
00937       
00938     using Teuchos::createVector;
00939 
00940     getGraph(br, bc, graphData, rowPtrs, nnzPerRow);
00941     ccmf->configure(lowestRow_[br], createVector(rowPtrs), createVector(nnzPerRow), createVector(graphData));
00942   }
00943   else
00944   {
00945     Tabs tab1;
00946     SUNDANCE_MSG2(verb, tab1 << "Assembler: doing incremental matrix structuring...");
00947     incrementalGetGraph(br, bc, icmf);
00948     {
00949       TimeMonitor timer1(matFinalizeTimer());
00950       icmf->finalize();
00951     }
00952   }
00953   
00954   SUNDANCE_MSG3(verb, tab << "Assembler: allocating matrix...");
00955   {
00956     TimeMonitor timer1(matAllocTimer());
00957     A = matFactory->createMatrix();
00958   }
00959 }
00960 
00961 Playa::LinearOperator<double> Assembler::allocateMatrix() const
00962 {
00963   LinearOperator<double> A;
00964   Array<Vector<double> > b;
00965   configureMatrix(A, b);
00966   return A;
00967 }
00968 
00969 
00970 
00971 
00972 
00973   
00974 void Assembler::displayEvaluationResults(
00975   const EvalContext& context, 
00976   const EvaluatableExpr* evalExpr, 
00977   const Array<double>& constantCoeffs, 
00978   const Array<RCP<EvalVector> >& vectorCoeffs) const 
00979 {
00980   Tabs tab;
00981   FancyOStream& os = Out::os();
00982 
00983   os << tab << "evaluation results: " << std::endl;
00984 
00985   const RCP<SparsitySuperset>& sparsity 
00986     = evalExpr->sparsitySuperset(context);
00987   
00988   sparsity->print(os, vectorCoeffs, constantCoeffs);
00989 }
00990 
00991 
00992 
00993 void Assembler::assemblyLoop(const ComputationType& compType,
00994   RCP<AssemblyKernelBase> kernel) const
00995 {
00996   Tabs tab;
00997   int verb = 0;
00998   if (eqn_->hasActiveWatchFlag()) verb = max(eqn_->maxWatchFlagSetting("assembly loop"), 1);
00999   
01000 
01001   SUNDANCE_BANNER1(verb, tab, "Assembly loop");
01002 
01003   SUNDANCE_MSG1(verb, tab << "computation type is " << compType); 
01004   /* Allocate space for the workset's list of cell local IDs.
01005    * Currently, a workset is an array of cell indices. It could be an array
01006    * of pointers to cell objects, cell iterators, or whatever is needed to
01007    * work with something like a Peano curve data structure. */
01008   SUNDANCE_MSG2(verb, tab << "work set size is " << workSetSize()); 
01009   RCP<Array<int> > workSet = rcp(new Array<int>());
01010   workSet->reserve(workSetSize());
01011 
01012   /* Allocate isLocalFlag array, which distinguishes between local
01013    * and non-local cells in the workset. This is used to prevent 
01014    * adding multiple copies of zero-form values for border cells. */
01015   RCP<Array<int> > isLocalFlag = rcp(new Array<int>());
01016   isLocalFlag->reserve(workSetSize());
01017 
01018   /* Declare objects for storage of coefficient evaluation results 
01019    * that are returned from the symbolic evaluation system. EvalVector
01020    * is the object in which a vector of numerical values for a 
01021    * given functional derivative are returned from an evaluation of
01022    * a symbolic DAG. */
01023   Array<RCP<EvalVector> > vectorCoeffs;
01024   Array<double> constantCoeffs;
01025 
01026   /* Create an object in which to store local integration results */
01027   RCP<Array<double> > localValues = rcp(new Array<double>());
01028 
01029   /* Get the symbolic specification of the current computation.
01030    * The "context" is simply a unique ID used to distinguish different
01031    * settings in which evaluation might be made. The same expression might be
01032    * evaluated in a context where some subset of all possible functional
01033    * derivatives are needed, and later in another context where a different 
01034    * subset of functional derivatives is needed. The "context" object lets
01035    * us associate a different Evaluator object with each such set of
01036    * requirements. */
01037   const Array<EvalContext>& contexts = contexts_.get(compType);
01038   /* Get the integral groups needed for this calculation */
01039   const Array<Array<RCP<IntegralGroup> > >& groups = groups_.get(compType);
01040   /* Get the expressions needed for this calculation */
01041   const Array<const EvaluatableExpr*>& evalExprs 
01042     = evalExprs_.get(compType);
01043   /* Get the flags indicating shich, if any, RQCs are to be skipped in this
01044    * calculation */
01045   const Array<int>& skipRqc = skipRqc_.get(compType);
01046 
01047   
01048   /* === Start main loop. 
01049    * The outer loop is over RQCs, that is, over unique combinations of subregions
01050    * (CellFilters) and quadrature rules.   */
01051   SUNDANCE_MSG1(verb, 
01052     tab << "---------- outer assembly loop over subregions");
01053   //SUNDANCE_MSG3(verb, tab << "Row DoF:" << rowMap_.size());
01054   //SUNDANCE_MSG3(verb, tab << "Column DoF:" << colMap_.size());
01055   //SUNDANCE_MSG3(verb, tab << "Region Quadratic Comb:" << rqc_.size());
01056 
01057   /* The double array which contains the transformations*/
01058   Array< Array<RCP<AssemblyTransformationBuilder> > > transformations;
01059 
01060   /** ----- Create Transformation objects for each integral group ------- */
01061   transformations.resize(rqc_.size());
01062   for (int r=0; r<rqc_.size(); r++)
01063   {
01064     transformations[r].resize(groups[r].size());
01065     for (int g=0; g<groups[r].size(); g++)
01066     {
01067       const RCP<IntegralGroup>& group = groups[r][g];
01068       /* Here we create the transformation object, if they are not needed
01069        * there would be no operation done to the array of local stiffness matrix */
01070       transformations[r][g] = rcp(new AssemblyTransformationBuilder( group , rowMap_ , colMap_ , mesh_));
01071     }
01072   }
01073 
01074 
01075   /* Record the default kernel verbosity so that it if changes we can
01076    * reset it at the end of a loop iteration */
01077   int oldKernelVerb = kernel->verb();
01078   
01079   TEUCHOS_TEST_FOR_EXCEPT(rqc_.size() != evalExprs.size());
01080 
01081   /* Looping over RQCs */
01082   for (int r=0; r<rqc_.size(); r++)
01083   {
01084     Tabs tab0;
01085 
01086     /* Set the verbosity level for this RQC */
01087     int rqcVerb = 0;
01088     int evalVerb = 0;
01089     int evalMedVerb = 0;
01090     int dfEvalVerb = 0;
01091     int fillVerb = 0;
01092 
01093     /* Check for watch point, and set verbosity accordingly */
01094     if (rqc_[r].watch().isActive()) 
01095     {
01096       Tabs tab01;
01097       rqcVerb=verb;
01098       evalVerb = rqc_[r].watch().param("evaluation");
01099       evalMedVerb = rqc_[r].watch().param("eval mediator");
01100       dfEvalVerb = rqc_[r].watch().param("discrete function evaluation");
01101       fillVerb = rqc_[r].watch().param("fill");
01102 
01103       SUNDANCE_MSG1(rqcVerb, tab0 << std::endl 
01104         << tab0 << "-------------"
01105         << std::endl << tab0 << " doing watched subregion r=" << r << " of " 
01106         << rqc_.size() << ", rqc=" 
01107         << rqc_[r]);    
01108       if (skipRqc[r]) 
01109       {
01110         SUNDANCE_MSG2(rqcVerb, tab01 << "this rqc is not needed in comp type = " << compType);
01111       }
01112       else
01113       {
01114         SUNDANCE_MSG2(rqcVerb, tab01 << "expr is " << evalExprs[r]->toString());
01115         SUNDANCE_MSG2(rqcVerb, tab01 << "isBC= " << isBCRqc_[r]); 
01116       }
01117     }
01118     else
01119     {
01120       SUNDANCE_MSG1(verb, tab0 << "unwatched region r=" << r << " of " << rqc_.size());
01121     }
01122     Tabs tab01;
01123 
01124     kernel->setVerb(fillVerb);
01125     
01126     /* Deciding whether we should skip this RQC in the current computation 
01127      * type. For example, a certain boundary surface might appear in the
01128      * computation of a functional but not in the state equation. */
01129     if ((!isBCRqc_[r] && eqn_->skipRqc(compType, rqc_[r]))
01130       || (isBCRqc_[r] && eqn_->skipBCRqc(compType, rqc_[r])))
01131     {
01132       Tabs tab012;
01133       SUNDANCE_MSG2(rqcVerb, tab012 << "nothing to do for comp type " 
01134         << compType);
01135       continue;
01136     }
01137 
01138     /* specify the evaluation mediator for this RQC.
01139      * Recall that the evaluation mediator is the object responsible for communication
01140      * between the symbolic expression tree and discretization-dependent data structures
01141      * such as discrete functions. 
01142      *
01143      * The evaluation manager is an object that is responsible for management
01144      * of the symbolic evaluation; it controls allocation of memory for evaluation
01145      * results, access to the evaluation mediator, and other administrative tasks. 
01146      */
01147     evalMgr_->setMediator(mediators_[r]);
01148     mediators_[r]->setVerb(evalMedVerb, dfEvalVerb);
01149 
01150     /* Tell the manager which CellFilter and QuadratureFamily we're currently working with. 
01151      * This is simply forwarded to the mediator, which needs to know the number
01152      * of quadrature points as well as mesh properties such as cell dimension. */
01153     evalMgr_->setRegion(contexts_.get(compType)[r]);
01154   
01155     /* get the cell filter for the current RQC */
01156     CellFilter filter = rqc_[r].domain();
01157     /* Find the cells that "pass" the predicate in the filter. Note: a CellFilter
01158      * will cache the cell sets it has computed, so the predicate computation will
01159      * only be done once per mesh, regardless of how often this function is called. 
01160      * Todo: the cache needs to be reset upon mesh refinement. 
01161      */
01162     CellSet cells = filter.getCells(mesh_);
01163     /* Find the dimension of cells that pass the current filter. */
01164     int cellDim = filter.dimension(mesh_);
01165     /* Find the type of cells in the current filter. Note: we've assumed here that all
01166      * cells have identical topology, and need to work out how to deal with 
01167      * meshes with mxed cell types. That will usually be handled by grouping similar
01168      * cells into "blocks" (as is done in Exodus files, for instance) in which case 
01169      * will still work, but there should be an error check to ensure that that assumption
01170      * is never violated. */ 
01171     CellType cellType = mesh_.cellType(cellDim);
01172     /* Get the cell type of the maximal-dimension cofacets, in case we need 
01173      * integrals or DOF maps done on the maximal cofacets. Todo: this code will break
01174      * for internal boundaries at the interface between cofacets of different types,
01175      * e.g., a face joining a prism and a hex. Not sure how to handle that case. */
01176     CellType maxCellType = mesh_.cellType(mesh_.spatialDim());
01177 
01178     SUNDANCE_MSG2(rqcVerb, tab01 << "cell type = " << cellType 
01179       << ", cellDim = " << cellDim 
01180       << ", max cell type = " << maxCellType 
01181       << ", max cell dim = " << mesh_.spatialDim());
01182 
01183 
01184     /* Determine whether we need to refer to maximal cofacets for 
01185      * some or all integrations and DOF mappings. */
01186     IntegrationCellSpecifier intCellSpec
01187       = rqcRequiresMaximalCofacets_.get(compType)[r];
01188     SUNDANCE_MSG2(rqcVerb, tab01 
01189       << "whether we need to refer to maximal cofacets: " << intCellSpec);
01190 
01191     /* Find the unknowns and variations appearing on the current domain. This
01192      * information is stored in the EquationSet object.  */
01193     const Array<Set<int> >& requiredVars = eqn_->reducedVarsOnRegion(filter);
01194     const Array<Set<int> >& requiredUnks = eqn_->reducedUnksOnRegion(filter);
01195 
01196     /* Prepare for evaluation on the current domain:
01197      * Tell the mediator whether maximal cofacets should be used (which will determine
01198      * how discrete functions are computed). */
01199     SUNDANCE_MSG2(rqcVerb, tab01 << "setting integration specifier in mediator");
01200     mediators_[r]->setIntegrationSpec(intCellSpec);
01201     /*
01202      * Tell the mediator the cell type, and whether we are on an internal
01203      * boundary. We need to know if we're on an internal boundary so that 
01204      * we can use DOF lookup logic that's capable of figuring out which
01205      * of two cofacets contains DOF information for those functions defined
01206      * on only one side of the boundary (as in, e.g., fluid-structure boundaries).
01207      */
01208     SUNDANCE_MSG2(rqcVerb, tab01 << "setting cell type=" << cellType << " in mediator");
01209     mediators_[r]->setCellType(cellType, maxCellType, isInternalBdry_[r]);    
01210     /* Get the Evaluator object that will actually carry out calculations on
01211      * the expression DAG in the current context (recall that a single expression
01212      * may support multiple evaluators). */
01213     const Evaluator* evaluator 
01214       = evalExprs[r]->evaluator(contexts[r]).get();
01215 
01216     /* Loop over cells in batches of the work set size.
01217      * At present, we're accumulating cell indices into an array. That would
01218      * need to be changed to work with Peano. */
01219     CellIterator iter=cells.begin();
01220     int workSetCounter = 0;
01221     int myRank = mesh_.comm().getRank();
01222 
01223     SUNDANCE_MSG2(rqcVerb, tab01 << "----- looping over worksets");
01224     while (iter != cells.end())
01225     {
01226       Tabs tab1;
01227       /* build up the work set: add cells until the work set size is 
01228        * reached or we run out of cells. To begin with, empty both the
01229        * workset array and the isLocalFlag array. Note that the reserve()
01230        * method has been called previously, so that as we append cells
01231        * to the array, no memory allocation is done (unless we run over 
01232        * the reserved size). */
01233       workSet->resize(0);
01234       isLocalFlag->resize(0);
01235       for (int c=0; c<workSetSize() && iter != cells.end(); c++, iter++)
01236       {
01237         workSet->append(*iter);
01238         /* we need the isLocalFlag values so that we can ignore contributions
01239          * to zero-forms from off-processor elements */
01240         isLocalFlag->append(myRank==mesh_.ownerProcID(cellDim, *iter));
01241       }
01242       /* The work set has now been accumulated */
01243       SUNDANCE_MSG2(rqcVerb,
01244         tab1 << "doing work set " << workSetCounter
01245         << " consisting of " << workSet->size() << " cells");
01246       {
01247         Tabs tab2;
01248         SUNDANCE_MSG4(rqcVerb, tab2 << "cells are " << *workSet);
01249       }
01250       workSetCounter++;
01251 
01252       /* set the verbosity for the evaluation mediator */
01253       evalMgr_->setVerb(evalVerb);
01254 
01255       /* Register the workset with the mediator. Internally, the mediator
01256        * will look up the cell Jacobians and facet indices needed for this calculation. It 
01257        * uses them for discrete function transformation. */
01258       mediators_[r]->setCellBatch(workSet);
01259       /* Get the Jacobians from the mediator, so that we can use the same Jacobian
01260        * objects for discrete function transformation and for integral 
01261        * transformation. The "volume" Jacobian is used to scale the integration
01262        * cell volume by det(J). The "transformation" Jacobian is used to
01263        * transform vectors. These will be the same, except for the case
01264        * where we integrate on a facet but do transformations by reference to 
01265        * a maximal cofacet. */
01266       const CellJacobianBatch& JVol = mediators_[r]->JVol();
01267       const CellJacobianBatch& JTrans = mediators_[r]->JTrans();
01268       /* Get from the mediator the facet indices for each cell. If I am integrating
01269        * on a facet (e.g., a boundary cell) but getting DOFs or JTrans from
01270        * a maximal cofacet, I need to know my index in the array of 
01271        * that cofacet's facets. */
01272       const Array<int>& facetIndices = mediators_[r]->facetIndices();
01273       if (facetIndices.size() > 0)
01274   {
01275     Tabs tab2;
01276     SUNDANCE_MSG2(rqcVerb, tab2 << "facet indices are " << facetIndices);
01277   }
01278 
01279       /* Reset the assembly kernel for the current workset. What happens at this
01280        * step depends on the specific kernel being used. The kernel might, for instance,
01281        * build local DOF maps for the current batch of cells. */
01282       kernel->prepareForWorkSet(
01283         requiredVars, 
01284         requiredUnks, 
01285         mediators_[r]);
01286         
01287       /* Evaluate the coefficient expressions. Recall that each coefficient
01288        * appearing in an integral is a particular functional derivative of the
01289        * integrand. The constant-valued coefficients are written into the
01290        * constantCoeffs array, the variable coefficients into the vectorCoeffs
01291        * array. 
01292        *
01293        * Recall that the eval manager contains the current evaluation mediator, so that
01294        * by passing the evaluation manager as an argument to evaluate(), the evaluation
01295        * is aware of the mediator and can therefore access discrete functions, etc.
01296        */ 
01297       evaluator->resetNumCalls();
01298       SUNDANCE_MSG2(rqcVerb, tab1 
01299         << "====== evaluating coefficient expressions") ;
01300       try
01301       {
01302         evalExprs[r]->evaluate(*evalMgr_, constantCoeffs, vectorCoeffs);
01303       }
01304       catch(std::exception& exc)
01305       {
01306         Tabs tabX;
01307         SUNDANCE_BANNER1(10, tabX, 
01308           "DETECTED EXCEPTION DURING EXPR EVALUATION CALL IN ASSEMBLY LOOP");
01309         Tabs tabX1;
01310         SUNDANCE_MSG1(10, tabX1 << "While working on RQC="
01311           << rqc_[r]);
01312         SUNDANCE_MSG1(10, tabX1 << "While evaluating expr="
01313           << evalExprs[r]->toString());
01314         throw (std::runtime_error(exc.what()));
01315       }
01316 
01317       /* Optionally, print the evaluation results */
01318       SUNDANCE_MSG2(rqcVerb, tab1 << "====== done evaluating expressions") ;
01319       if (evalVerb > 3) 
01320       {
01321         displayEvaluationResults(contexts[r], evalExprs[r], constantCoeffs, 
01322           vectorCoeffs);
01323       }
01324     
01325       /* ---- Do the element integrals and insertions ------ */
01326 
01327       /* The matrices used to transform integrals are built upon first use by 
01328        * this workset, then cached because they may be needed for several 
01329        * integrals on the same workset. As we're now in a new workset with
01330        * new cells, they should be rebuilt if needed. This step informs all integrals
01331        * that the cache is out of date. 
01332        *
01333        * Todo: this uses a static function that contains static data (via the "Meyers
01334        * trick") and so is not thread-safe. If we want to do multithreaded parallelism
01335        * for multicore architectures, this implementation will need to be changed.
01336        */ 
01337       ElementIntegral::invalidateTransformationMatrices();
01338       
01339       /* Loop over the integral groups */
01340       SUNDANCE_MSG2(rqcVerb, tab1 << "----- looping over integral groups");
01341       for (int g=0; g<groups[r].size(); g++)
01342       {
01343         Tabs tab2;
01344         SUNDANCE_MSG2(rqcVerb, tab2 << std::endl << tab2 
01345           << "--- evaluating integral group g=" << g << " of " 
01346           << groups[r].size() );
01347 
01348         /* Do the integrals. The integration results will be written into
01349          * the array "localValues". */
01350         const RCP<IntegralGroup>& group = groups[r][g];
01351         if (!group->evaluate(JTrans, JVol, *isLocalFlag, facetIndices, workSet,
01352             vectorCoeffs, constantCoeffs, localValues)) continue;
01353 
01354         /* Here we call the transformation object, if they are not needed
01355          * (the function might be one return) there would be no operation
01356          * done to the array of local stiffness matrix
01357          * Do the actual transformation (transformations for Matrix)*/
01358         transformations[r][g]->applyTransformsToAssembly(
01359             g , (localValues->size() / workSet->size()),
01360                 cellType, cellDim , maxCellType,
01361             JTrans, JVol, facetIndices, workSet, localValues);
01362 
01363         /* add the integration results into the output objects by a call
01364          * to the kernel's fill() function. We need to pass isBCRqc to the kernel
01365          * because it might handle BC rows differently. The integral group
01366          * data structure contains information about which test and unknown
01367          * functions were used in this integration, and so provides to the assembly
01368          * kernel such information as is needed to look up the correct DOFs for this
01369          * batch of integrals. */
01370         {
01371           TimeMonitor fillTM(fillTimer());
01372           kernel->fill(isBCRqc_[r], *group, localValues);
01373         }
01374       }
01375       SUNDANCE_MSG2(rqcVerb, tab1 << "----- done looping over integral groups");
01376     }
01377     SUNDANCE_MSG2(rqcVerb, tab0 << "----- done looping over worksets");
01378     /* reset the kernel verbosity to the default */
01379     kernel->setVerb(oldKernelVerb);
01380     SUNDANCE_MSG1(verb, tab0 << "----- done rqc");
01381   }
01382   SUNDANCE_MSG1(verb, tab << "----- done looping over rqcs");
01383 
01384 
01385   /* Do any post-fill processing, such as MPI_AllReduce add on functional values. */
01386   SUNDANCE_MSG2(verb, tab << "doing post-loop processing"); 
01387   kernel->postLoopFinalization();
01388   SUNDANCE_BANNER1(verb, tab, "done assembly loop"); 
01389 
01390   /* All done with assembly! */
01391 }
01392 
01393 
01394 
01395 /* ------------  assemble both the vector and the matrix  ------------- */
01396 
01397 void Assembler::assemble(LinearOperator<double>& A,
01398   Array<Vector<double> >& mv) const 
01399 {
01400   TimeMonitor timer(assemblyTimer());
01401   Tabs tab;
01402   int verb = 0;
01403   if (eqn_->hasActiveWatchFlag()) verb = max(verb, 1);
01404   
01405   SUNDANCE_BANNER1(verb, tab, "Assembling matrix and vector");
01406 
01407   TEUCHOS_TEST_FOR_EXCEPTION(!contexts_.containsKey(MatrixAndVector),
01408     std::runtime_error,
01409     "Assembler::assemble(A, b) called for an assembler that "
01410     "does not support matrix/vector assembly");
01411 
01412   configureMatrix(A, mv);
01413   
01414   RCP<AssemblyKernelBase> kernel 
01415     = rcp(new MatrixVectorAssemblyKernel(
01416             rowMap_, isBCRow_, lowestRow_,
01417             colMap_, isBCCol_, lowestCol_,
01418             A, mv, partitionBCs_, 
01419             0));
01420 
01421   assemblyLoop(MatrixAndVector, kernel);
01422 
01423   SUNDANCE_MSG1(verb, tab << "matrix=" << A);
01424   if (verb>0) A.print(Out::os());
01425   SUNDANCE_MSG1(verb, tab << "vectors=" << mv);
01426   for (int i=0; i<mv.size(); i++) 
01427   {
01428     SUNDANCE_MSG1(verb, tab << "vectors #" << i);
01429     if (verb>0) mv[i].print(Out::os());
01430   }
01431 
01432   SUNDANCE_MSG1(verb, tab << "Assembler: done assembling matrix and vector");
01433 }
01434 
01435 /* ------------  assemble the matrix and sensitivity RHS ------------- */
01436 
01437 void Assembler::assembleSensitivities(LinearOperator<double>& A,
01438   Array<Vector<double> >& mv) const 
01439 {
01440   TimeMonitor timer(assemblyTimer());
01441   Tabs tab;
01442   int verb = 0;
01443   if (eqn_->hasActiveWatchFlag()) verb = max(verb, 1);
01444   
01445   SUNDANCE_BANNER1(verb, tab, "Assembling matrix and sensitivity vector");
01446 
01447   TEUCHOS_TEST_FOR_EXCEPTION(!contexts_.containsKey(Sensitivities),
01448     std::runtime_error,
01449     "Assembler::assembleSensitivities(A, b) called for an assembler that "
01450     "does not support sensitivity assembly");
01451 
01452   configureMatrix(A, mv);
01453   
01454   
01455   RCP<AssemblyKernelBase> kernel 
01456     = rcp(new MatrixVectorAssemblyKernel(
01457             rowMap_, isBCRow_, lowestRow_,
01458             colMap_, isBCCol_, lowestCol_,
01459             A, mv, partitionBCs_, 
01460             0));
01461 
01462   assemblyLoop(Sensitivities, kernel);
01463   SUNDANCE_MSG1(verb, tab << "Assembler: done assembling matrix and sensitivity vector");
01464 }
01465 
01466 
01467 /* ------------  assemble the vector alone  ------------- */
01468 
01469 void Assembler::assemble(Array<Vector<double> >& mv) const 
01470 {
01471   /* Tab is advanced by ctor, taken back by dtor upon leaving scope */
01472   Tabs tab;
01473   /* Timer is started by ctor, stopped by dtor upon leaving scope */
01474   TimeMonitor timer(assemblyTimer());  
01475 
01476   /* If any subexpression is watched, print basic information */ 
01477   int verb = 0;
01478   if (eqn_->hasActiveWatchFlag()) verb = max(verb, 1);
01479 
01480   SUNDANCE_BANNER1(verb, tab, "Assembling vector");
01481 
01482   /* Throw an exception if we don't know how to compute a vector */
01483   TEUCHOS_TEST_FOR_EXCEPTION(!contexts_.containsKey(VectorOnly),
01484     std::runtime_error,
01485     "Assembler::assemble(b) called for an assembler that "
01486     "does not support vector-only assembly");
01487 
01488   /* The vector is not configured yet. Do so here. */
01489   configureVector(mv);
01490   
01491   /* Create an assembly kernel that knows how to fill a vector */
01492   RCP<AssemblyKernelBase> kernel 
01493     = rcp(new VectorAssemblyKernel(
01494             rowMap_, isBCRow_, lowestRow_,
01495             mv, partitionBCs_, 0));
01496 
01497   assemblyLoop(VectorOnly, kernel);
01498 
01499   SUNDANCE_MSG1(verb, tab << "Assembler: done assembling vector");
01500 }
01501 
01502 
01503 /* ------------  evaluate a functional and its gradient ---- */
01504 
01505 void Assembler::evaluate(double& value, Array<Vector<double> >& gradient) const 
01506 {
01507   Tabs tab;
01508   TimeMonitor timer(assemblyTimer());
01509   int verb = 0;
01510   if (eqn_->hasActiveWatchFlag()) verb = max(verb, 1);
01511 
01512   SUNDANCE_BANNER1(verb, tab, "Computing functional and gradient");
01513 
01514   TEUCHOS_TEST_FOR_EXCEPTION(!contexts_.containsKey(FunctionalAndGradient),
01515     std::runtime_error,
01516     "Assembler::evaluate(f,df) called for an assembler that "
01517     "does not support value/gradient assembly");
01518 
01519   configureVector(gradient);
01520 
01521   value = 0.0;
01522   
01523   RCP<AssemblyKernelBase> kernel 
01524     = rcp(new FunctionalGradientAssemblyKernel(
01525             mesh_.comm(),
01526             rowMap_, isBCRow_, lowestRow_,
01527             gradient, partitionBCs_, &value, verb));
01528 
01529   assemblyLoop(FunctionalAndGradient, kernel);
01530 
01531   SUNDANCE_BANNER1(verb, tab, "Done computing functional and gradient");
01532 
01533 }
01534 
01535 
01536 
01537 
01538 /* ------------  evaluate a functional ---- */
01539 
01540 void Assembler::evaluate(double& value) const 
01541 {
01542   Tabs tab;
01543   TimeMonitor timer(assemblyTimer());
01544   int verb = 0;
01545   if (eqn_->hasActiveWatchFlag()) verb = max(verb, 1);
01546 
01547   SUNDANCE_BANNER1(verb, tab, "Computing functional");
01548 
01549   TEUCHOS_TEST_FOR_EXCEPTION(!contexts_.containsKey(FunctionalOnly),
01550     std::runtime_error,
01551     "Assembler::evaluate(f) called for an assembler that "
01552     "does not support functional evaluation");
01553 
01554   value = 0.0;
01555   
01556   RCP<AssemblyKernelBase> kernel 
01557     = rcp(new FunctionalAssemblyKernel(mesh_.comm(), &value, verb));
01558 
01559   assemblyLoop(FunctionalOnly, kernel);
01560 
01561   SUNDANCE_BANNER1(verb, tab, "Done computing functional");
01562 }
01563 
01564 
01565 
01566 
01567 /* ------------  get the nonzero pattern for the matrix ------------- */
01568                        
01569                        
01570 void Assembler::getGraph(int br, int bc, 
01571   Array<int>& graphData,
01572   Array<int>& rowPtrs,
01573   Array<int>& nnzPerRow) const 
01574 {
01575   TimeMonitor timer(graphBuildTimer());
01576   Tabs tab;
01577   int verb = eqn_->maxWatchFlagSetting("matrix config");
01578 
01579 
01580   RCP<Array<int> > workSet = rcp(new Array<int>());
01581   workSet->reserve(workSetSize());
01582 
01583   RCP<Array<Array<int> > > testLocalDOFs 
01584     = rcp(new Array<Array<int> >());
01585 
01586   RCP<Array<Array<int> > > unkLocalDOFs
01587     = rcp(new Array<Array<int> >());
01588 
01589   SUNDANCE_MSG3(verb, tab << "Creating graph: there are " 
01590     << rowMap_[br]->numLocalDOFs()
01591     << " local equations");
01592 
01593 
01594   Array<Set<int> > tmpGraph;
01595   tmpGraph.resize(rowMap_[br]->numLocalDOFs());
01596 
01597   {
01598     TimeMonitor timer2(colSearchTimer());
01599     for (int d=0; d<eqn_->numRegions(); d++)
01600     {
01601       Tabs tab0;
01602       CellFilter domain = eqn_->region(d);
01603       const Array<Set<int> >& requiredVars = eqn_->reducedVarsOnRegion(domain);
01604       const Array<Set<int> >& requiredUnks = eqn_->reducedUnksOnRegion(domain);
01605       SUNDANCE_MSG2(verb,
01606         tab0 << "cell set " << domain
01607         << " isBCRegion=" << eqn_->isBCRegion(d));
01608       int dim = domain.dimension(mesh_);
01609       CellSet cells = domain.getCells(mesh_);
01610 
01611       RCP<Set<OrderedPair<int, int> > > pairs ;
01612       if (eqn_->hasVarUnkPairs(domain)) pairs = eqn_->varUnkPairs(domain);
01613 
01614       SUNDANCE_OUT(verb > 2 && pairs.get() != 0, 
01615         tab0 << "non-BC pairs = "
01616         << *pairs);
01617        
01618       RCP<Set<OrderedPair<int, int> > > bcPairs ;
01619       if (eqn_->isBCRegion(d))
01620       {
01621         if (eqn_->hasBCVarUnkPairs(domain)) 
01622         {
01623           bcPairs = eqn_->bcVarUnkPairs(domain);
01624           SUNDANCE_OUT(verb > 2, tab0 << "BC pairs = "
01625             << *bcPairs);
01626         }
01627       }
01628       Array<Set<int> > unksForTestsSet(eqn_->numVarIDs(bc));
01629       Array<Set<int> > bcUnksForTestsSet(eqn_->numVarIDs(bc));
01630 
01631       Set<OrderedPair<int, int> >::const_iterator i;
01632       
01633       if (pairs.get() != 0)
01634       {
01635         for (i=pairs->begin(); i!=pairs->end(); i++)
01636         {
01637           const OrderedPair<int, int>& p = *i;
01638           int t = p.first();
01639           int u = p.second();
01640 
01641           TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasVarID(t), std::logic_error,
01642             "Test function ID " << t << " does not appear "
01643             "in equation set");
01644           TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasUnkID(u), std::logic_error,
01645             "Unk function ID " << u << " does not appear "
01646             "in equation set");
01647 
01648           if (eqn_->blockForVarID(t) != br) continue;
01649           if (eqn_->blockForUnkID(u) != bc) continue;
01650 
01651           unksForTestsSet[eqn_->reducedVarID(t)].put(eqn_->reducedUnkID(u));
01652         }
01653       }
01654       if (bcPairs.get() != 0)
01655       {
01656         for (i=bcPairs->begin(); i!=bcPairs->end(); i++)
01657         {
01658           const OrderedPair<int, int>& p = *i;
01659           int t = p.first();
01660           int u = p.second();
01661 
01662           if (eqn_->blockForVarID(t) != br) continue;
01663           if (eqn_->blockForUnkID(u) != bc) continue;
01664 
01665           TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasVarID(t), std::logic_error,
01666             "Test function ID " << t << " does not appear "
01667             "in equation set");
01668           TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasUnkID(u), std::logic_error,
01669             "Unk function ID " << u << " does not appear "
01670             "in equation set");
01671           bcUnksForTestsSet[eqn_->reducedVarID(t)].put(eqn_->reducedUnkID(u));
01672         }
01673       }
01674 
01675       Array<Array<int> > unksForTests(unksForTestsSet.size());
01676       Array<Array<int> > bcUnksForTests(bcUnksForTestsSet.size());
01677 
01678       for (int t=0; t<unksForTests.size(); t++)
01679       {
01680         unksForTests[t] = unksForTestsSet[t].elements();
01681         bcUnksForTests[t] = bcUnksForTestsSet[t].elements();
01682       }
01683       
01684       Array<int> numTestNodes;
01685       Array<int> numUnkNodes;
01686 
01687       
01688 
01689       int highestRow = lowestRow_[br] + rowMap_[br]->numLocalDOFs();
01690 
01691       int nt = eqn_->numVarIDs(br);
01692 
01693       CellIterator iter=cells.begin();
01694       while (iter != cells.end())
01695       {
01696         /* build a work set */
01697         workSet->resize(0);
01698         for (int c=0; c<workSetSize() && iter != cells.end(); c++, iter++)
01699         {
01700           workSet->append(*iter);
01701         }
01702 
01703         int nCells = workSet->size();
01704 
01705         RCP<const MapStructure> colMapStruct; 
01706         RCP<const MapStructure> rowMapStruct 
01707           = rowMap_[br]->getDOFsForCellBatch(dim, *workSet, 
01708             requiredVars[br], *testLocalDOFs,
01709             numTestNodes, verb);
01710         if (rowMap_[br].get()==colMap_[bc].get())
01711         {
01712           unkLocalDOFs = testLocalDOFs;
01713           numUnkNodes = numTestNodes;
01714           colMapStruct = rowMapStruct;
01715         }
01716         else
01717         {
01718           colMapStruct = colMap_[br]->getDOFsForCellBatch(dim, *workSet, 
01719             requiredUnks[bc], 
01720             *unkLocalDOFs, numUnkNodes, verb);
01721         }
01722 
01723         if (pairs.get() != 0)
01724         {
01725           for (int c=0; c<nCells; c++)
01726           {
01727             for (int t=0; t<nt; t++)
01728             {
01729               int tChunk = rowMapStruct->chunkForFuncID(t);
01730               int nTestFuncs = rowMapStruct->numFuncs(tChunk);
01731               int testFuncIndex = rowMapStruct->indexForFuncID(t);
01732               int nTestNodes = numTestNodes[tChunk];
01733               const Array<int>& testDOFs = (*testLocalDOFs)[tChunk];
01734               for (int uit=0; uit<unksForTests[t].size(); uit++)
01735               {
01736                 Tabs tab2;
01737                 int u = unksForTests[t][uit];
01738                 int uChunk = colMapStruct->chunkForFuncID(u);
01739                 int nUnkFuncs = colMapStruct->numFuncs(uChunk);
01740                 int unkFuncIndex = colMapStruct->indexForFuncID(u);
01741                 const Array<int>& unkDOFs = (*unkLocalDOFs)[uChunk];
01742                 int nUnkNodes = numUnkNodes[uChunk];
01743                 for (int n=0; n<nTestNodes; n++)
01744                 {
01745                   int row
01746                     = testDOFs[(c*nTestFuncs + testFuncIndex)*nTestNodes + n];
01747                   if (row < lowestRow_[br] || row >= highestRow
01748                     || (*(isBCRow_[br]))[row-lowestRow_[br]]) continue;
01749                   Set<int>& colSet = tmpGraph[row-lowestRow_[br]];
01750                   for (int m=0; m<nUnkNodes; m++)
01751                   {
01752                     int col 
01753                       = unkDOFs[(c*nUnkFuncs + unkFuncIndex)*nUnkNodes + m];
01754                     colSet.put(col);
01755                   }
01756                 }
01757               }
01758             }
01759           }
01760         }
01761         if (bcPairs.get() != 0)
01762         {
01763           for (int c=0; c<nCells; c++)
01764           {
01765             for (int t=0; t<nt; t++)
01766             {
01767               int tChunk = rowMapStruct->chunkForFuncID(t);
01768               int nTestFuncs = rowMapStruct->numFuncs(tChunk);
01769               int testFuncIndex = rowMapStruct->indexForFuncID(t);
01770               int nTestNodes = numTestNodes[tChunk];
01771               const Array<int>& testDOFs = (*testLocalDOFs)[tChunk];
01772               for (int uit=0; uit<bcUnksForTests[t].size(); uit++)
01773               {
01774                 Tabs tab2;
01775                 int u = bcUnksForTests[t][uit];
01776                 int uChunk = colMapStruct->chunkForFuncID(u);
01777                 int nUnkFuncs = colMapStruct->numFuncs(uChunk);
01778                 int unkFuncIndex = colMapStruct->indexForFuncID(u);
01779                 const Array<int>& unkDOFs = (*unkLocalDOFs)[uChunk];
01780                 int nUnkNodes = numUnkNodes[uChunk];
01781                 for (int n=0; n<nTestNodes; n++)
01782                 {
01783                   int row
01784                     = testDOFs[(c*nTestFuncs + testFuncIndex)*nTestNodes + n];
01785                   if (row < lowestRow_[br] || row >= highestRow
01786                     || !(*(isBCRow_[br]))[row-lowestRow_[br]]) continue;
01787                   Set<int>& colSet = tmpGraph[row-lowestRow_[br]];
01788                   for (int m=0; m<nUnkNodes; m++)
01789                   {
01790                     int col 
01791                       = unkDOFs[(c*nUnkFuncs + unkFuncIndex)*nUnkNodes + m];
01792                     colSet.put(col);
01793                   }
01794                 }
01795               }
01796             }
01797           }
01798         }
01799       }
01800     }
01801   }
01802 
01803   
01804   {
01805     TimeMonitor t2(graphFlatteningTimer());
01806     int nLocalRows = rowMap_[br]->numLocalDOFs();
01807 
01808     int nnz = 0;
01809     rowPtrs.resize(nLocalRows);
01810     nnzPerRow.resize(rowMap_[br]->numLocalDOFs());
01811     for (int i=0; i<nLocalRows; i++) 
01812     {
01813       rowPtrs[i] = nnz;
01814       nnzPerRow[i] = tmpGraph[i].size();
01815       nnz += nnzPerRow[i];
01816     }
01817 
01818     graphData.resize(nnz);
01819     int* base = &(graphData[0]);
01820     for (int i=0; i<nLocalRows; i++)
01821     {
01822       //        tmpGraph[i].fillArray(base + rowPtrs[i]);
01823       int* rowBase = base + rowPtrs[i];
01824       const Set<int>& rowSet = tmpGraph[i];
01825       int k = 0;
01826       for (Set<int>::const_iterator 
01827              j=rowSet.begin(); j != rowSet.end(); j++, k++)
01828       {
01829         rowBase[k] = *j;
01830       }
01831     }
01832   }
01833 
01834 }
01835 /* ------------  get the nonzero pattern for the matrix ------------- */
01836                        
01837                        
01838 void Assembler
01839 ::incrementalGetGraph(int br, int bc,
01840   IncrementallyConfigurableMatrixFactory* icmf) const 
01841 {
01842   TimeMonitor timer(graphBuildTimer());
01843   Tabs tab;
01844   int verb = eqn_->maxWatchFlagSetting("matrix config");
01845 
01846   RCP<Array<int> > workSet = rcp(new Array<int>());
01847   workSet->reserve(workSetSize());
01848 
01849   RCP<Array<Array<int> > > testLocalDOFs 
01850     = rcp(new Array<Array<int> >());
01851 
01852   RCP<Array<Array<int> > > unkLocalDOFs
01853     = rcp(new Array<Array<int> >());
01854 
01855   SUNDANCE_MSG2(verb, tab << "Creating graph: there are " 
01856     << rowMap_[br]->numLocalDOFs()
01857     << " local equations");
01858 
01859 
01860   for (int d=0; d<eqn_->numRegions(); d++)
01861   {
01862     Tabs tab0;
01863     CellFilter domain = eqn_->region(d);
01864     const Array<Set<int> >& requiredVars = eqn_->reducedVarsOnRegion(domain);
01865     const Array<Set<int> >& requiredUnks = eqn_->reducedUnksOnRegion(domain);
01866     Array<int> localVars = requiredVars[br].elements();
01867     Array<int> localUnks = requiredUnks[bc].elements();
01868     SUNDANCE_MSG3(verb,
01869       tab0 << "cell set " << domain
01870       << " isBCRegion=" << eqn_->isBCRegion(d));
01871     int dim = domain.dimension(mesh_);
01872     CellSet cells = domain.getCells(mesh_);
01873 
01874     RCP<Set<OrderedPair<int, int> > > pairs ;
01875     if (eqn_->hasVarUnkPairs(domain)) pairs = eqn_->varUnkPairs(domain);
01876 
01877     SUNDANCE_OUT(verb > 2 && pairs.get() != 0, 
01878       tab0 << "non-BC pairs = "
01879       << *pairs);
01880        
01881     RCP<Set<OrderedPair<int, int> > > bcPairs ;
01882     if (eqn_->isBCRegion(d))
01883     {
01884       if (eqn_->hasBCVarUnkPairs(domain)) 
01885       {
01886         bcPairs = eqn_->bcVarUnkPairs(domain);
01887         SUNDANCE_MSG3(verb, tab0 << "BC pairs = "
01888           << *bcPairs);
01889       }
01890     }
01891     Array<Set<int> > unksForTestsSet(eqn_->numVarIDs(br));
01892     Array<Set<int> > bcUnksForTestsSet(eqn_->numVarIDs(br));
01893 
01894     Set<OrderedPair<int, int> >::const_iterator i;
01895       
01896     if (pairs.get() != 0)
01897     {
01898       for (i=pairs->begin(); i!=pairs->end(); i++)
01899       {
01900         const OrderedPair<int, int>& p = *i;
01901         int t = p.first();
01902         int u = p.second();
01903 
01904         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasVarID(t), std::logic_error,
01905           "Test function ID " << t << " does not appear "
01906           "in equation set");
01907         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasUnkID(u), std::logic_error,
01908           "Unk function ID " << u << " does not appear "
01909           "in equation set");
01910 
01911 
01912         if (eqn_->blockForVarID(t) != br) continue;
01913         if (eqn_->blockForUnkID(u) != bc) continue;
01914 
01915         unksForTestsSet[eqn_->reducedVarID(t)].put(eqn_->reducedUnkID(u));
01916       }
01917     }
01918     if (bcPairs.get() != 0)
01919     {
01920       for (i=bcPairs->begin(); i!=bcPairs->end(); i++)
01921       {
01922         const OrderedPair<int, int>& p = *i;
01923         int t = p.first();
01924         int u = p.second();
01925         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasVarID(t), std::logic_error,
01926           "Test function ID " << t << " does not appear "
01927           "in equation set");
01928         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasUnkID(u), std::logic_error,
01929           "Unk function ID " << u << " does not appear "
01930           "in equation set");
01931 
01932         if (eqn_->blockForVarID(t) != br) continue;
01933         if (eqn_->blockForUnkID(u) != bc) continue;
01934 
01935         bcUnksForTestsSet[eqn_->reducedVarID(t)].put(eqn_->reducedUnkID(u));
01936       }
01937     }
01938 
01939     Array<Array<int> > unksForTests(unksForTestsSet.size());
01940     Array<Array<int> > bcUnksForTests(bcUnksForTestsSet.size());
01941 
01942     for (int t=0; t<unksForTests.size(); t++)
01943     {
01944       unksForTests[t] = unksForTestsSet[t].elements();
01945       bcUnksForTests[t] = bcUnksForTestsSet[t].elements();
01946     }
01947       
01948     Array<int> numTestNodes;
01949     Array<int> numUnkNodes;
01950       
01951     int highestRow = lowestRow_[br] + rowMap_[br]->numLocalDOFs();
01952 
01953     CellIterator iter=cells.begin();
01954 
01955     while (iter != cells.end())
01956     {
01957       /* build a work set */
01958       workSet->resize(0);
01959       for (int c=0; c<workSetSize() && iter != cells.end(); c++, iter++)
01960       {
01961         workSet->append(*iter);
01962       }
01963 
01964       int nCells = workSet->size();
01965 
01966       RCP<const MapStructure> colMapStruct; 
01967       RCP<const MapStructure> rowMapStruct 
01968         = rowMap_[br]->getDOFsForCellBatch(dim, *workSet, 
01969           requiredVars[br],
01970           *testLocalDOFs,
01971           numTestNodes, verb);
01972 
01973       if (rowMap_[br].get()==colMap_[bc].get())
01974       {
01975         unkLocalDOFs = testLocalDOFs;
01976         numUnkNodes = numTestNodes;
01977         colMapStruct = rowMapStruct;
01978       }
01979       else
01980       {
01981         colMapStruct = colMap_[bc]->getDOFsForCellBatch(dim, *workSet, 
01982           requiredUnks[bc],
01983           *unkLocalDOFs, numUnkNodes, verb);
01984       }
01985 
01986           
01987       if (pairs.get() != 0)
01988       {
01989         for (int c=0; c<nCells; c++)
01990         {
01991           for (int tIndex=0; tIndex<localVars.size(); tIndex++)
01992           {
01993             int t = localVars[tIndex];
01994             int tChunk = rowMapStruct->chunkForFuncID(t);
01995             int nTestFuncs = rowMapStruct->numFuncs(tChunk);
01996             int testFuncIndex = rowMapStruct->indexForFuncID(t);
01997             int nTestNodes = numTestNodes[tChunk];
01998             const Array<int>& testDOFs = (*testLocalDOFs)[tChunk];
01999             for (int uit=0; uit<unksForTests[t].size(); uit++)
02000             {
02001               Tabs tab2;
02002               int u = unksForTests[t][uit];
02003               int uChunk = colMapStruct->chunkForFuncID(u);
02004               int nUnkFuncs = colMapStruct->numFuncs(uChunk);
02005               int unkFuncIndex = colMapStruct->indexForFuncID(u);
02006               const Array<int>& unkDOFs = (*unkLocalDOFs)[uChunk];
02007               int nUnkNodes = numUnkNodes[uChunk];
02008               for (int n=0; n<nTestNodes; n++)
02009               {
02010                 int row
02011                   = testDOFs[(c*nTestFuncs + testFuncIndex)*nTestNodes + n];
02012                 if (row < lowestRow_[br] || row >= highestRow
02013                   || (*(isBCRow_[br]))[row-lowestRow_[br]]) continue;
02014                 const int* colPtr = &(unkDOFs[(c*nUnkFuncs + unkFuncIndex)*nUnkNodes]);
02015                 icmf->initializeNonzerosInRow(row, nUnkNodes, colPtr);
02016               }
02017             }
02018           }
02019         }
02020       }
02021       if (bcPairs.get() != 0)
02022       {
02023         for (int c=0; c<nCells; c++)
02024         {
02025           for (int tIndex=0; tIndex<localVars.size(); tIndex++)
02026           {
02027             int t = localVars[tIndex];
02028             int tChunk = rowMapStruct->chunkForFuncID(t);
02029             int nTestFuncs = rowMapStruct->numFuncs(tChunk);
02030             int testFuncIndex = rowMapStruct->indexForFuncID(t);
02031             int nTestNodes = numTestNodes[tChunk];
02032             const Array<int>& testDOFs = (*testLocalDOFs)[tChunk];
02033             for (int uit=0; uit<bcUnksForTests[t].size(); uit++)
02034             {
02035               Tabs tab2;
02036               int u = bcUnksForTests[t][uit];
02037               int uChunk = colMapStruct->chunkForFuncID(u);
02038               int nUnkFuncs = colMapStruct->numFuncs(uChunk);
02039               int unkFuncIndex = colMapStruct->indexForFuncID(u);
02040               const Array<int>& unkDOFs = (*unkLocalDOFs)[uChunk];
02041               int nUnkNodes = numUnkNodes[uChunk];
02042               for (int n=0; n<nTestNodes; n++)
02043               {
02044                 int row
02045                   = testDOFs[(c*nTestFuncs + testFuncIndex)*nTestNodes + n];
02046                 if (row < lowestRow_[br] || row >= highestRow
02047                   || !(*(isBCRow_[br]))[row-lowestRow_[br]]) continue;
02048 
02049                 const int* colPtr = &(unkDOFs[(c*nUnkFuncs + unkFuncIndex)*nUnkNodes]);
02050                 icmf->initializeNonzerosInRow(row, nUnkNodes, colPtr);
02051               }
02052             }
02053           }
02054         }
02055       }
02056     }
02057   }
02058 }
02059 
02060 
02061 Array<Array<int> > Assembler::findNonzeroBlocks() const
02062 {
02063   int verb = eqn_->maxWatchFlagSetting("assembler setup");
02064 
02065   Array<Array<int> > rtn(eqn_->numVarBlocks());
02066   for (int br=0; br<rtn.size(); br++)
02067   {
02068     rtn[br].resize(eqn_->numUnkBlocks());
02069     for (int bc=0; bc<rtn[br].size(); bc++)
02070     {
02071       rtn[br][bc] = 0 ;
02072     }
02073   }
02074 
02075   for (int d=0; d<eqn_->numRegions(); d++)
02076   {
02077     Tabs tab0;
02078     CellFilter domain = eqn_->region(d);
02079     SUNDANCE_MSG3(verb,
02080       tab0 << "cell set " << domain
02081       << " isBCRegion=" << eqn_->isBCRegion(d));
02082 
02083     RCP<Set<OrderedPair<int, int> > > pairs ;
02084     if (eqn_->hasVarUnkPairs(domain)) pairs = eqn_->varUnkPairs(domain);
02085 
02086     SUNDANCE_OUT(verb > 2 && pairs.get() != 0, 
02087       tab0 << "non-BC pairs = "
02088       << *pairs);
02089        
02090     RCP<Set<OrderedPair<int, int> > > bcPairs ;
02091     if (eqn_->isBCRegion(d))
02092     {
02093       if (eqn_->hasBCVarUnkPairs(domain)) 
02094       {
02095         bcPairs = eqn_->bcVarUnkPairs(domain);
02096         SUNDANCE_MSG3(verb, tab0 << "BC pairs = "
02097           << *bcPairs);
02098       }
02099     }
02100 
02101     Set<OrderedPair<int, int> >::const_iterator i;
02102       
02103     if (pairs.get() != 0)
02104     {
02105       for (i=pairs->begin(); i!=pairs->end(); i++)
02106       {
02107         const OrderedPair<int, int>& p = *i;
02108         int t = p.first();
02109         int u = p.second();
02110 
02111         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasVarID(t), std::logic_error,
02112           "Test function ID " << t << " does not appear "
02113           "in equation set");
02114         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasUnkID(u), std::logic_error,
02115           "Unk function ID " << u << " does not appear "
02116           "in equation set");
02117 
02118 
02119         int br = eqn_->blockForVarID(t);
02120         int bc = eqn_->blockForUnkID(u);
02121         rtn[br][bc] = 1;
02122       }
02123     }
02124     if (bcPairs.get() != 0)
02125     {
02126       for (i=bcPairs->begin(); i!=bcPairs->end(); i++)
02127       {
02128         const OrderedPair<int, int>& p = *i;
02129         int t = p.first();
02130         int u = p.second();
02131         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasVarID(t), std::logic_error,
02132           "Test function ID " << t << " does not appear "
02133           "in equation set");
02134         TEUCHOS_TEST_FOR_EXCEPTION(!eqn_->hasUnkID(u), std::logic_error,
02135           "Unk function ID " << u << " does not appear "
02136           "in equation set");
02137         int br = eqn_->blockForVarID(t);
02138         int bc = eqn_->blockForUnkID(u);
02139         rtn[br][bc] = 1;
02140       }
02141     }
02142   }
02143 
02144   return rtn;
02145 }
02146 
02147 VectorSpace<double> Assembler::solnVecSpace() const
02148 {
02149   Array<VectorSpace<double> > rtn(eqn_->numUnkBlocks());
02150 
02151   for (int i=0; i<rtn.size(); i++)
02152   {
02153     rtn[i] = solutionSpace()[i]->vecSpace();
02154   }
02155 
02156   if ((int) rtn.size() == 1)
02157   {
02158     return rtn[0];
02159   }
02160   return blockSpace(rtn);
02161 }
02162 
02163 
02164 VectorSpace<double> Assembler::rowVecSpace() const
02165 {
02166   Array<VectorSpace<double> > rtn(eqn_->numVarBlocks());
02167 
02168   for (int i=0; i<rtn.size(); i++)
02169   {
02170     rtn[i] = rowSpace()[i]->vecSpace();
02171   }
02172 
02173   if ((int) rtn.size() == 1)
02174   {
02175     return rtn[0];
02176   }
02177   return blockSpace(rtn);
02178 }
02179 
02180 
02181 
02182 int& Assembler::workSetSize()
02183 {
02184   static int rtn = defaultWorkSetSize(); 
02185   return rtn;
02186 }
02187 
02188 int Assembler::maxWatchFlagSetting(const std::string& name) const 
02189 {
02190   return eqnSet()->maxWatchFlagSetting(name);
02191 }
02192 
02193 bool Assembler::matNeedsConfiguration() const
02194 {
02195   return Teuchos::is_null(cachedAssembledMatrix_.ptr());
02196 }
02197 
02198 
02199 void Assembler::flushConfiguration() const
02200 {
02201     numConfiguredColumns_ = 0;
02202     cachedAssembledMatrix_ = LinearOperator<double>();
02203     mesh_.flushSpecialWeights();
02204     mesh_.flushCurvePoints();
02205 
02206     // empty the cache from the cell filters
02207     for (int r=0; r<rqc_.size(); r++)
02208     {
02209       const CellFilter filter = rqc_[r].domain();
02210       filter.cfbPtr()->flushCache();
02211     }
02212 }

Site Contact