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

Site Contact