SundanceQuadratureEvalMediator.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 "SundanceQuadratureEvalMediator.hpp"
00032 #include "SundanceCoordExpr.hpp"
00033 #include "SundanceTempStack.hpp"
00034 #include "SundanceCellDiameterExpr.hpp"
00035 #include "SundanceCellVectorExpr.hpp"
00036 #include "SundanceSpatialDerivSpecifier.hpp"
00037 #include "SundanceDiscreteFunction.hpp"
00038 #include "SundanceDiscreteFuncElement.hpp"
00039 #include "SundanceCellJacobianBatch.hpp"
00040 #include "SundanceOut.hpp"
00041 #include "SundanceTabs.hpp"
00042 #include "SundanceExceptions.hpp"
00043 
00044 #include "Teuchos_BLAS.hpp"
00045 
00046 
00047 using namespace Sundance;
00048 using namespace Sundance;
00049 using namespace Sundance;
00050 using namespace Sundance;
00051 using namespace Sundance;
00052 using namespace Sundance;
00053 using namespace Sundance;
00054 using namespace Teuchos;
00055 using namespace TSFExtended;
00056 
00057 TEUCHOS_TIMER(coordEvalTimer, "Quad mediator: coord eval")
00058 
00059 
00060 QuadratureEvalMediator
00061 ::QuadratureEvalMediator(const Mesh& mesh, 
00062   int cellDim,
00063   const QuadratureFamily& quad)
00064   : StdFwkEvalMediator(mesh, cellDim),
00065     numEvaluationCases_(-1),
00066     quad_(quad),
00067     numQuadPtsForCellType_(),
00068     quadPtsForReferenceCell_(),
00069     quadPtsReferredToMaxCell_(),
00070     physQuadPts_(),
00071     refFacetBasisVals_(2)
00072 {}
00073 
00074 Time& QuadratureEvalMediator::coordEvaluationTimer()
00075 {
00076   return coordEvalTimer();
00077 }
00078 
00079 void QuadratureEvalMediator::setCellType(const CellType& cellType,
00080   const CellType& maxCellType,
00081   bool isInternBdry) 
00082 {
00083   StdFwkEvalMediator::setCellType(cellType, maxCellType, isInternBdry);
00084 
00085   Tabs tab;
00086 
00087   SUNDANCE_MSG2(verb(), tab <<  "QuadEvalMed::setCellType: cellType=" 
00088     << cellType << " cellDim=" << cellDim() << " maxCellType=" << maxCellType);
00089   SUNDANCE_MSG2(verb(), tab << "integration spec: =" 
00090     << integrationCellSpec());
00091   SUNDANCE_MSG2(verb(), tab << "forbid cofacet integrations =" 
00092     << forbidCofacetIntegrations());
00093   if (isInternalBdry()) 
00094   {
00095     SUNDANCE_MSG2(verb(), tab << "working on internal boundary");
00096   }
00097   
00098 //  TEST_FOR_EXCEPT(isInternalBdry()
00099 //    && integrationCellSpec() != NoTermsNeedCofacets);
00100 
00101   if (cellType != maxCellType)
00102   {
00103     Tabs tab1;
00104     SUNDANCE_MSG2(verb(), tab1 << "working out #facet cases"); 
00105     numEvaluationCases_ = numFacets(maxCellType, cellDim());
00106   }
00107   else 
00108   {
00109     Tabs tab1;
00110     SUNDANCE_MSG2(verb(), tab1 << "no need for facet cases; work on original cell"); 
00111     numEvaluationCases_ = 1;
00112   }
00113   SUNDANCE_MSG2(verb(), tab << "num eval cases =" << numEvaluationCases_); 
00114 
00115   if (!isInternalBdry() 
00116     && quadPtsReferredToMaxCell_.containsKey(cellType)) return;
00117 
00118   if (!quadPtsForReferenceCell_.containsKey(cellType))
00119   {
00120     SUNDANCE_MSG2(verb(), tab << "creating quad points for ref cell type=" 
00121       << cellType);
00122     RCP<Array<Point> > pts = rcp(new Array<Point>());
00123     RCP<Array<double> > wgts = rcp(new Array<double>()); 
00124     
00125     quad_.getPoints(cellType, *pts, *wgts);
00126     quadPtsForReferenceCell_.put(cellType, pts);
00127     
00128     numQuadPtsForCellType_.put(cellType, pts->size());
00129   }
00130 
00131   if (!quadPtsReferredToMaxCell_.containsKey(cellType))
00132   {
00133     SUNDANCE_MSG2(verb(), tab << "creating quad points for max cell type=" 
00134       << maxCellType);
00135     RCP<Array<Array<Point> > > facetPts 
00136       = rcp(new Array<Array<Point> >(numEvaluationCases()));
00137     RCP<Array<Array<double> > > facetWgts 
00138       = rcp(new Array<Array<double> >(numEvaluationCases()));
00139 
00140     for (int fc=0; fc<numEvaluationCases(); fc++)
00141     {
00142       if (cellType != maxCellType)
00143       {
00144         quad_.getFacetPoints(maxCellType, cellDim(), fc, 
00145           (*facetPts)[fc], (*facetWgts)[fc]);
00146       }
00147       else
00148       {
00149         quad_.getPoints(maxCellType, (*facetPts)[fc], (*facetWgts)[fc]);
00150       }
00151     }
00152     quadPtsReferredToMaxCell_.put(cellType, facetPts);
00153   }
00154 }
00155 
00156 int QuadratureEvalMediator::numQuadPts(const CellType& ct) const 
00157 {
00158   TEST_FOR_EXCEPTION(!numQuadPtsForCellType_.containsKey(ct),
00159     RuntimeError,
00160     "no quadrature points have been tabulated for cell type=" << ct);
00161   return numQuadPtsForCellType_.get(ct);
00162 }
00163 
00164 void QuadratureEvalMediator::evalCellDiameterExpr(const CellDiameterExpr* expr,
00165   RCP<EvalVector>& vec) const 
00166 {
00167   Tabs tabs;
00168   SUNDANCE_MSG2(verb(),tabs 
00169     << "QuadratureEvalMediator evaluating cell diameter expr " 
00170     << expr->toString());
00171 
00172   int nQuad = numQuadPts(cellType());
00173   int nCells = cellLID()->size();
00174 
00175   SUNDANCE_MSG3(verb(),tabs << "number of quad pts=" << nQuad);
00176   Array<double> diameters;
00177   mesh().getCellDiameters(cellDim(), *cellLID(), diameters);
00178 
00179   vec->resize(nQuad*nCells);
00180   double * const xx = vec->start();
00181   int k=0;
00182   for (int c=0; c<nCells; c++)
00183   {
00184     double h = diameters[c];
00185     for (int q=0; q<nQuad; q++, k++) 
00186     {
00187       xx[k] = h;
00188     }
00189   }
00190 }
00191 
00192 void QuadratureEvalMediator::evalCellVectorExpr(const CellVectorExpr* expr,
00193   RCP<EvalVector>& vec) const 
00194 {
00195   Tabs tabs;
00196   SUNDANCE_MSG2(verb(),tabs 
00197     << "QuadratureEvalMediator evaluating cell vector expr " 
00198     << expr->toString());
00199 
00200   int nQuad = numQuadPts(cellType());
00201   int nCells = cellLID()->size();
00202 
00203   vec->resize(nQuad*nCells);
00204 
00205   SUNDANCE_MSG3(verb(),tabs << "number of quad pts=" << nQuad);
00206   int dir = expr->componentIndex();
00207 
00208   Array<Point> vectors;
00209   if (expr->isNormal())
00210   { 
00211     mesh().outwardNormals(*cellLID(), vectors);
00212   }
00213   else
00214   {
00215     TEST_FOR_EXCEPTION(cellDim() != 1, RuntimeError,
00216       "unable to compute tangent vectors for cell dim = " << cellDim());
00217     mesh().tangentsToEdges(*cellLID(), vectors);
00218   }
00219     
00220   double * const xx = vec->start();
00221   int k=0;
00222   for (int c=0; c<nCells; c++)
00223   {
00224     double n = vectors[c][dir];
00225     for (int q=0; q<nQuad; q++, k++) 
00226     {
00227       xx[k] = n;
00228     }
00229   }
00230 }
00231 
00232 void QuadratureEvalMediator::evalCoordExpr(const CoordExpr* expr,
00233   RCP<EvalVector>& vec) const 
00234 {
00235   Tabs tabs;
00236   SUNDANCE_MSG2(verb(),tabs 
00237     << "QuadratureEvalMediator evaluating coord expr " 
00238     << expr->toString());
00239 
00240   TimeMonitor timer(coordEvalTimer());
00241   
00242   computePhysQuadPts();
00243   int nQuad = physQuadPts_.length();
00244   int d = expr->dir();
00245   
00246   SUNDANCE_MSG3(verb(),tabs << "number of quad pts=" << nQuad);
00247 
00248   vec->resize(nQuad);
00249   double * const xx = vec->start();
00250   for (int q=0; q<nQuad; q++) 
00251   {
00252     xx[q] = physQuadPts_[q][d];
00253   }
00254 }
00255 
00256 RCP<Array<Array<Array<double> > > > QuadratureEvalMediator
00257 ::getFacetRefBasisVals(const BasisFamily& basis, int diffOrder) const
00258 {
00259   Tabs tab;
00260   RCP<Array<Array<Array<double> > > > rtn ;
00261 
00262   CellType evalCellType = cellType();
00263   if (cellDim() != maxCellDim())
00264   {
00265     if (verb() >= 2)
00266     {
00267       Out::os() << tab << "alwaysUseCofacets = " 
00268                 << ElementIntegral::alwaysUseCofacets() << std::endl;
00269       Out::os() << tab << "diffOrder = " << diffOrder << std::endl;
00270     }
00271     if (ElementIntegral::alwaysUseCofacets() || diffOrder>0)
00272     {
00273       if (!cofacetCellsAreReady()) setupFacetTransformations();
00274       evalCellType = maxCellType();
00275     
00276       TEST_FOR_EXCEPTION(!cofacetCellsAreReady(), RuntimeError, 
00277         "cofacet cells not ready in getFacetRefBasisVals()");
00278     }
00279   }
00280 
00281   SUNDANCE_MSG2(verb(), tab << "eval cell type = " << evalCellType);
00282   typedef OrderedPair<BasisFamily, CellType> key;
00283 
00284   int nDerivResults = 1;
00285   if (diffOrder==1) nDerivResults = maxCellDim();
00286 
00287   if (!refFacetBasisVals_[diffOrder].containsKey(key(basis, cellType())))
00288   {
00289     SUNDANCE_OUT(this->verb() > 2,
00290       tab << "computing basis values on facet quad pts");
00291     rtn = rcp(new Array<Array<Array<double> > >(numEvaluationCases()));
00292 
00293     Array<Array<Array<Array<double> > > > tmp(nDerivResults);
00294     
00295     if (verb() >= 2)
00296     {
00297       Out::os() << tab << "numEvalCases = " << numEvaluationCases()
00298                 << std::endl;
00299       Out::os() << tab << "diff order = " << diffOrder << std::endl;
00300       Out::os() << tab << "cell type = " << cellType() << std::endl;
00301       Out::os() << tab << "quad pt map = ";
00302       if (evalCellType!=cellType())
00303       { 
00304         Out::os() << quadPtsReferredToMaxCell_ << std::endl;
00305       }
00306       else
00307       {
00308         Out::os() << quadPtsForReferenceCell_ << std::endl;
00309       }
00310     }
00311 
00312     if (evalCellType!=cellType())
00313     { 
00314       TEST_FOR_EXCEPTION(quadPtsReferredToMaxCell_.size() == 0,
00315         RuntimeError,
00316         "empty quadrature point map (max cell)");
00317     }
00318     else
00319     {
00320       TEST_FOR_EXCEPTION(quadPtsForReferenceCell_.size() == 0,
00321         RuntimeError,
00322         "empty quadrature point map (ref cell)");
00323     }
00324 
00325     for (int fc=0; fc<numEvaluationCases(); fc++)
00326     {
00327       Tabs tab1;
00328       (*rtn)[fc].resize(basis.dim());
00329       SUNDANCE_MSG2(verb(), tab1 << "fc = " << fc);
00330 
00331       for (int r=0; r<nDerivResults; r++)
00332       {
00333         tmp[r].resize(basis.dim());
00334         MultiIndex mi;
00335         if (diffOrder==1)
00336         {
00337           mi[r]=1;
00338         }
00339         SpatialDerivSpecifier deriv(mi);
00340         /* Here we evaluate the basis functions at specified quadrature points
00341          * on the reference cell */
00342         if (evalCellType != cellType())
00343         {
00344           SUNDANCE_MSG2(verb(), tab1 << "referring to max cell");
00345           basis.refEval(evalCellType, 
00346             (*(quadPtsReferredToMaxCell_.get(cellType())))[fc], 
00347             deriv, tmp[r], verb());
00348         }
00349         else
00350         {
00351           SUNDANCE_MSG2(verb(), tab1 << "computing on reference cell");
00352           basis.refEval(evalCellType, 
00353             (*(quadPtsForReferenceCell_.get(cellType()))), 
00354             deriv, tmp[r], verb());
00355         }
00356       }
00357       /* the tmp array contains values indexed as [quad][node]. 
00358        * We need to put this into fortran order with quad index running
00359        * fastest */
00360       int dim = maxCellDim();
00361       int nQuad = tmp[0][0].size();
00362       int nNodes = tmp[0][0][0].size();
00363       int nTot = dim * nQuad * nNodes;
00364       for (int d=0; d<basis.dim(); d++)
00365       {
00366         (*rtn)[fc][d].resize(nTot);
00367         for (int r=0; r<nDerivResults; r++)
00368         {
00369           for (int q=0; q<nQuad; q++)
00370           {
00371             for (int n=0; n<nNodes; n++)
00372             {
00373               (*rtn)[fc][d][(n*nQuad + q)*nDerivResults + r] = tmp[r][d][q][n];
00374             }
00375           }
00376         }
00377       }
00378     }
00379     refFacetBasisVals_[diffOrder].put(key(basis, cellType()), rtn);
00380   }
00381   else
00382   {
00383     SUNDANCE_OUT(this->verb() > 2,
00384       tab << "reusing facet basis values on quad pts");
00385     rtn = refFacetBasisVals_[diffOrder].get(key(basis, cellType()));
00386   }
00387 
00388   return rtn;
00389 }
00390 
00391 Array<Array<double> >* QuadratureEvalMediator
00392 ::getRefBasisVals(const BasisFamily& basis, int diffOrder) const
00393 {
00394   Tabs tab;
00395   RCP<Array<Array<Array<double> > > > fRtn 
00396     = getFacetRefBasisVals(basis, diffOrder);
00397   return &((*fRtn)[0]);
00398 }
00399 
00400 
00401 void QuadratureEvalMediator
00402 ::evalDiscreteFuncElement(const DiscreteFuncElement* expr,
00403   const Array<MultiIndex>& multiIndices,
00404   Array<RCP<EvalVector> >& vec) const
00405 {
00406   const DiscreteFunctionData* f = DiscreteFunctionData::getData(expr);
00407   TEST_FOR_EXCEPTION(f==0, InternalError,
00408     "QuadratureEvalMediator::evalDiscreteFuncElement() called "
00409     "with expr that is not a discrete function");
00410   Tabs tab;
00411 
00412   SUNDANCE_MSG1(verb(),tab << "QuadEvalMed evaluating DF " << expr->name());
00413 
00414   int nQuad = numQuadPts(cellType());
00415   int myIndex = expr->myIndex();
00416 
00417   for (int i=0; i<multiIndices.size(); i++)
00418   {
00419     Tabs tab1;
00420     const MultiIndex& mi = multiIndices[i];
00421     SUNDANCE_MSG2(dfVerb(),
00422       tab1 << "evaluating DF " << expr->name() 
00423       << " for multiindex " << mi << std::endl
00424       << tab1 << "num cells = " << cellLID()->size() << std::endl
00425       << tab1 << "num quad points = " << nQuad << std::endl
00426       << tab1 << "my index = " << expr->myIndex() << std::endl
00427       << tab1 << "num funcs = " << f->discreteSpace().nFunc());
00428 
00429     vec[i]->resize(cellLID()->size() * nQuad);
00430   
00431     if (mi.order() == 0)
00432     {
00433       Tabs tab2;
00434       SUNDANCE_MSG2(dfVerb(),tab2 << "in mi.order() == 0 branch");
00435       if (!fCache().containsKey(f) || !fCacheIsValid()[f])
00436       {
00437         fillFunctionCache(f, mi);
00438       }
00439       else
00440       {
00441         SUNDANCE_MSG2(dfVerb(),tab2 << "reusing function cache");
00442       }
00443 
00444       const RCP<const MapStructure>& mapStruct = mapStructCache()[f];
00445       int chunk = mapStruct->chunkForFuncID(myIndex);
00446       int funcIndex = mapStruct->indexForFuncID(myIndex);
00447       int nFuncs = mapStruct->numFuncs(chunk);
00448 
00449       SUNDANCE_MSG3(dfVerb(),tab2 << "chunk number = " << chunk << std::endl
00450         << tab2 << "function index=" << funcIndex << " of nFuncs=" 
00451         << nFuncs);
00452 
00453       const RCP<Array<Array<double> > >& cacheVals 
00454         = fCache()[f];
00455 
00456       SUNDANCE_MSG4(dfVerb(),tab2 << "cached function values=" << (*cacheVals)[chunk]);
00457 
00458       const double* cachePtr = &((*cacheVals)[chunk][0]);
00459       double* vecPtr = vec[i]->start();
00460           
00461       int cellSize = nQuad*nFuncs;
00462       int offset = funcIndex*nQuad;
00463       SUNDANCE_MSG3(dfVerb(),tab2 << "cell size=" << cellSize << ", offset=" 
00464         << offset);
00465       int k = 0;
00466       for (int c=0; c<cellLID()->size(); c++)
00467       {
00468         for (int q=0; q<nQuad; q++, k++)
00469         {
00470           vecPtr[k] = cachePtr[c*cellSize + offset + q];
00471         }
00472       }
00473       SUNDANCE_MSG4(dfVerb(),tab2 << "result vector=");
00474       if (dfVerb() >= 5)
00475       {
00476         vec[i]->print(Out::os());
00477         computePhysQuadPts();
00478         k=0;
00479         for (int c=0; c<cellLID()->size(); c++)
00480         {
00481           Out::os() << tab2 << "-------------------------------------------"
00482                     << std::endl;
00483           Out::os() << tab2 << "c=" << c << " cell LID=" << (*cellLID())[c]
00484                     << std::endl;
00485           Tabs tab3;
00486           for (int q=0; q<nQuad; q++, k++)
00487           {
00488             Out::os() << tab3 << "q=" << q << " " << physQuadPts_[k]
00489                       << " val=" << vecPtr[k] << std::endl;
00490           }
00491         }
00492       }
00493     }
00494     else
00495     {
00496       Tabs tab2;
00497       SUNDANCE_MSG2(dfVerb(),tab2 << "in mi.order() != 0 branch");
00498       if (!dfCache().containsKey(f) || !dfCacheIsValid()[f])
00499       {
00500         fillFunctionCache(f, mi);
00501       }
00502       else
00503       {
00504         SUNDANCE_MSG3(dfVerb(),tab2 << "reusing function cache");
00505       }
00506 
00507       RCP<const MapStructure> mapStruct = mapStructCache()[f];
00508       int chunk = mapStruct->chunkForFuncID(myIndex);
00509       int funcIndex = mapStruct->indexForFuncID(myIndex);
00510       int nFuncs = mapStruct->numFuncs(chunk);
00511 
00512 
00513       SUNDANCE_MSG3(dfVerb(),tab2 << "chunk number = " << chunk << std::endl
00514         << tab2 << "function index=" << funcIndex << " of nFuncs=" 
00515         << nFuncs);
00516 
00517       const RCP<Array<Array<double> > >& cacheVals 
00518         = dfCache()[f];
00519 
00520       SUNDANCE_MSG4(dfVerb(),tab2 << "cached function values=" << (*cacheVals)[chunk]);
00521 
00522       int dim = maxCellDim();
00523       int pDir = mi.firstOrderDirection();
00524       const double* cachePtr = &((*cacheVals)[chunk][0]);
00525       double* vecPtr = vec[i]->start();
00526 
00527       int cellSize = nQuad*nFuncs*dim;
00528       int offset = funcIndex * nQuad * dim;
00529       int k = 0;
00530 
00531       SUNDANCE_MSG2(dfVerb(),tab2 << "dim=" << dim << ", pDir=" << pDir
00532         << ", cell size=" << cellSize << ", offset=" 
00533         << offset);
00534       for (int c=0; c<cellLID()->size(); c++)
00535       {
00536         for (int q=0; q<nQuad; q++, k++)
00537         {
00538           vecPtr[k] = cachePtr[c*cellSize + offset + q*dim + pDir];
00539         }
00540       }
00541       SUNDANCE_MSG4(dfVerb(),tab2 << "result vector=");
00542       if (dfVerb() >= 5)
00543       {
00544         vec[i]->print(Out::os());
00545         computePhysQuadPts();
00546         k=0;
00547         for (int c=0; c<cellLID()->size(); c++)
00548         {
00549           Out::os() << tab2 << "-------------------------------------------"
00550                     << std::endl;
00551           Out::os() << tab2 << "c=" << c << " cell LID=" << (*cellLID())[c]
00552                     << std::endl;
00553           Tabs tab3;
00554           for (int q=0; q<nQuad; q++, k++)
00555           {
00556             Out::os() << tab3 << "q=" << q << " " << physQuadPts_[k]
00557                       << " val=" << vecPtr[k] << std::endl;
00558           }
00559         }
00560       }
00561     }
00562   }
00563 }
00564 
00565 void QuadratureEvalMediator::fillFunctionCache(const DiscreteFunctionData* f,
00566   const MultiIndex& mi) const 
00567 {
00568   Tabs tab0;
00569   
00570   SUNDANCE_MSG2(dfVerb(), tab0 << "QuadratureEvalMediator::fillFunctionCache()");
00571   SUNDANCE_MSG2(dfVerb(), tab0 << "multiIndex=" << mi);
00572   
00573   
00574   int diffOrder = mi.order();
00575   CellType evalCellType = cellType();
00576 
00577   int flops = 0;
00578   double jFlops = CellJacobianBatch::totalFlops();
00579 
00580   RCP<Array<Array<double> > > localValues;
00581   RCP<const MapStructure> mapStruct;
00582 
00583   Teuchos::BLAS<int,double> blas;
00584 
00585   {
00586     Tabs tab1;
00587     SUNDANCE_MSG2(dfVerb(), tab1 << "packing local values");
00588     if (!localValueCacheIsValid().containsKey(f) 
00589       || !localValueCacheIsValid().get(f))
00590     {
00591       Tabs tab2;
00592       SUNDANCE_MSG2(dfVerb(), tab2 << "filling cache");
00593       localValues = rcp(new Array<Array<double> >());
00594       if (cellDim() != maxCellDim())
00595       {
00596         if (dfVerb() >= 2)
00597         {
00598           Out::os() << tab2 << "alwaysUseCofacets = " 
00599                     << ElementIntegral::alwaysUseCofacets() << std::endl;
00600           Out::os() << tab2 << "diffOrder = " << diffOrder << std::endl;
00601         }
00602         if (ElementIntegral::alwaysUseCofacets() || diffOrder>0)
00603         {
00604           if (!cofacetCellsAreReady()) setupFacetTransformations();
00605           mapStruct = f->getLocalValues(maxCellDim(), *cofacetCellLID(), 
00606             *localValues);
00607           evalCellType = maxCellType();
00608         }
00609         else
00610         {
00611           mapStruct = f->getLocalValues(cellDim(), *cellLID(), 
00612             *localValues);
00613         }
00614       }
00615       else
00616       {
00617         mapStruct = f->getLocalValues(maxCellDim(), *cellLID(), *localValues);
00618       }
00619 
00620       TEST_FOR_EXCEPT(mapStruct.get() == 0);
00621       localValueCache().put(f, localValues);
00622       mapStructCache().put(f, mapStruct);
00623       localValueCacheIsValid().put(f, true);
00624     }
00625     else
00626     {
00627       Tabs tab2;
00628       SUNDANCE_MSG2(dfVerb(), tab2 << "reusing cache");
00629       localValues = localValueCache().get(f);
00630       mapStruct = mapStructCache().get(f);
00631     }
00632   }
00633 
00634   RCP<Array<Array<double> > > cacheVals;
00635 
00636   if (mi.order()==0)
00637   {
00638     if (fCache().containsKey(f))
00639     {
00640       cacheVals = fCache().get(f);
00641     }
00642     else
00643     {
00644       cacheVals = rcp(new Array<Array<double> >(mapStruct->numBasisChunks()));
00645       fCache().put(f, cacheVals);
00646     }
00647     fCacheIsValid().put(f, true);
00648   }
00649   else
00650   {
00651     if (dfCache().containsKey(f))
00652     {
00653       cacheVals = dfCache().get(f);
00654     }
00655     else
00656     {
00657       cacheVals = rcp(new Array<Array<double> >(mapStruct->numBasisChunks()));
00658       dfCache().put(f, cacheVals);
00659     }
00660     dfCacheIsValid().put(f, true);
00661   }
00662 
00663 
00664   
00665   for (int chunk=0; chunk<mapStruct->numBasisChunks(); chunk++)
00666   {
00667     Tabs tab1;
00668     SUNDANCE_MSG2(dfVerb(), tab1 << "processing dof map chunk=" << chunk
00669       << " of " << mapStruct->numBasisChunks());
00670     BasisFamily basis = rcp_dynamic_cast<BasisFamilyBase>(mapStruct->basis(chunk));
00671     SUNDANCE_MSG4(dfVerb(), tab1 << "basis=" << basis);
00672 
00673     int nFuncs = mapStruct->numFuncs(chunk);
00674     SUNDANCE_MSG2(dfVerb(), tab1 << "num funcs in this chunk=" << nFuncs);
00675     
00676 
00677     Array<double>& cache = (*cacheVals)[chunk];
00678 
00679     int nQuad = numQuadPts(cellType());
00680     int nCells = cellLID()->size();
00681     SUNDANCE_MSG2(dfVerb(), tab1 << "num quad points=" << nQuad);
00682     SUNDANCE_MSG2(dfVerb(), tab1 << "num cells=" << nCells);
00683 
00684     int nDir;
00685 
00686     if (mi.order()==1)
00687     {
00688       nDir = maxCellDim();
00689     }
00690     else
00691     {
00692       nDir = 1;
00693     }
00694     cache.resize(cellLID()->size() * nQuad * nDir * nFuncs);
00695 
00696       
00697     /* 
00698      * Sum over nodal values, which we can do with a matrix-matrix multiply
00699      * between the ref basis values and the local function values.
00700      *
00701      * There are two cases: (1) When we are evaluating 
00702      * on a facet, we must use different sets of reference function values
00703      * on the different facets. We must therefore loop over the evaluation
00704      * cells, using a vector of reference values chosen according to the
00705      * facet number of the current cell.  Let A be the
00706      * (nQuad*nDir)-by-(nNode) matrix of reference basis values for the
00707      * current cell's facet index and B be the (nNode)-by-(nFuncs) matrix of
00708      * function coefficient values for the current cell. Then C = A * B is
00709      * the (nQuad*nDir)-by-(nFunc) matrix of the function's derivative
00710      * values at the quadrature points in the current cell.  Each
00711      * matrix-matrix multiplication is done with a call to dgemm.
00712      *
00713      * (2) In other cases, we're either evaluating spatial derivatives on a
00714      * maximal cell or evaluating 0-order derivatives on a submaximal
00715      * cell. In these cases, all cells in the workset have the same
00716      * reference values. This lets us reuse the same matrix A on all matrix
00717      * multiplications, so that we can assemble one big
00718      * (nNode)-by-(nFuncs*nCells) matrix B and do all cells with a single
00719      * dgemm call to multiply A*B. The result C is then a single
00720      * (nQuad*nDir)-by-(nFuncs*nCells) matrix.
00721 
00722     */
00723     if (cellType() != evalCellType)
00724     {
00725       Tabs tab2;
00726       SUNDANCE_MSG2(dfVerb(), 
00727         tab2 << "evaluating by reference to maximal cell");
00728       
00729       RCP<Array<Array<Array<double> > > > refFacetBasisValues 
00730         = getFacetRefBasisVals(basis, mi.order());
00731       /* Note: even though we're ultimately not evaluating on 
00732        * maxCellType() here, use maxCellType() for both arguments
00733        * to nReferenceDOFs() because derivatives need to be
00734        * evaluated using all DOFs from the maximal cell, not just
00735        * those on the facet.
00736        */
00737       int nNodes = basis.nReferenceDOFsWithFacets(maxCellType(), maxCellType());
00738       int nRowsA = nQuad*nDir;
00739       int nColsA = nNodes;
00740       int nColsB = nFuncs; 
00741       int lda = nRowsA;
00742       int ldb = nNodes;
00743       int ldc = lda;
00744       double alpha = 1.0;
00745       double beta = 0.0;
00746 
00747       SUNDANCE_MSG2(dfVerb(), tab2 << "building transformation matrices cell-by-cell");
00748       int vecComp = 0;
00749       for (int c=0; c<nCells; c++)
00750       {
00751         int facetIndex = facetIndices()[c];
00752         double* A = &((*refFacetBasisValues)[facetIndex][vecComp][0]);
00753         double* B = &((*localValues)[chunk][c*nNodes*nFuncs]);
00754         double* C = &((*cacheVals)[chunk][c*nRowsA*nColsB]);
00755         blas.GEMM( Teuchos::NO_TRANS, Teuchos::NO_TRANS, nRowsA, nColsB, nColsA,
00756           alpha, A, lda, B, ldb, beta, C, ldc);
00757       }
00758     }
00759     else /* cellType() == evalCellType */
00760     {
00761       /* 
00762        * Sum over nodal values, which we can do with a matrix-matrix multiply
00763        * between the ref basis values and the local function values.
00764        */
00765       Tabs tab2;
00766       SUNDANCE_MSG2(dfVerb(), tab2 << "building batch transformation matrix");
00767 
00768       Array<Array<double> >* refBasisValues 
00769         = getRefBasisVals(basis, diffOrder);
00770       int nNodes = basis.nReferenceDOFsWithFacets(maxCellType(), cellType());
00771       int nRowsA = nQuad*nDir;
00772       int nColsA = nNodes;
00773       int nColsB = nFuncs*nCells; 
00774       int lda = nRowsA;
00775       int ldb = nNodes;
00776       int ldc = lda;
00777       double alpha = 1.0;
00778       double beta = 0.0;
00779       int vecComp = 0;
00780       if (dfVerb() >= 5)
00781       {
00782         Tabs tab3;
00783         Out::os() << tab2 << "Printing values at nodes" << std::endl;
00784         for (int c=0; c<nCells; c++)
00785         {
00786           Out::os() << tab3 << "-------------------------------------------"
00787                     << std::endl;
00788           Out::os() << tab3 << "c=" << c << " cell LID=" << (*cellLID())[c]
00789                     << std::endl;
00790           Tabs tab4;
00791           int offset = c*nNodes*nFuncs;
00792           for (int n=0; n<nNodes; n++)
00793           {
00794             Out::os() << tab4 << "n=" << n;
00795             for (int f=0; f<nFuncs; f++)
00796             {
00797               Out::os() << " " << (*localValues)[chunk][offset + n*nFuncs + f];
00798             }
00799             Out::os() << std::endl;
00800           }
00801         }
00802         Out::os() << tab2 << "-------------------------------------------";
00803         Out::os() << tab2 << "Printing reference basis at nodes" << std::endl;
00804         Out::os() << tab2 << "-------------------------------------------"
00805                   << std::endl;
00806         for (int n=0; n<nNodes; n++)
00807         {
00808           Out::os() << tab3 << "node=" << n << std::endl;
00809           for (int q=0; q<nQuad; q++)
00810           {
00811             Tabs tab4;
00812             Out::os() << tab4 << "q=" << q;
00813             for (int r=0; r<nDir; r++)
00814             {
00815               Out::os () << " " 
00816                          << (*refBasisValues)[vecComp][(n*nQuad + q)*nDir + r];
00817             }
00818             Out::os() << std::endl;
00819           }
00820         }
00821       }
00822       double* A = &((*refBasisValues)[vecComp][0]);
00823       double* B = &((*localValues)[chunk][0]);
00824       double* C = &((*cacheVals)[chunk][0]);
00825       blas.GEMM( Teuchos::NO_TRANS, Teuchos::NO_TRANS, nRowsA, nColsB, nColsA, alpha,
00826         A, lda, B, ldb, beta, C, ldc );
00827     }
00828 
00829 
00830     /* Transform derivatives to physical coordinates */
00831     const CellJacobianBatch& J = JTrans();
00832     double* C = &((*cacheVals)[chunk][0]);
00833     if (mi.order()==1)
00834     {
00835       SUNDANCE_MSG2(dfVerb(), tab1 << "doing transformations via dgemm");
00836       Tabs tab2;
00837       SUNDANCE_MSG2(dfVerb(), tab2 << "Jacobian batch nCells=" << J.numCells());
00838       SUNDANCE_MSG2(dfVerb(), tab2 << "Jacobian batch cell dim=" << J.cellDim());
00839       SUNDANCE_MSG2(dfVerb(), tab2 << "Jacobian batch spatial dim=" << J.spatialDim());
00840     
00841       int nRhs = nQuad * nFuncs;
00842       for (int c=0; c<cellLID()->size(); c++)
00843       {
00844         double* rhsPtr = &(C[(nRhs * nDir)*c]);
00845         J.applyInvJ(c, 0, rhsPtr, nRhs, false);
00846       }
00847     }
00848     else
00849     {
00850       SUNDANCE_MSG2(dfVerb(), tab1 << "no derivatives to transform");
00851     }
00852     SUNDANCE_MSG2(dfVerb(), tab1 << "done transformations");
00853   }
00854 
00855   jFlops = CellJacobianBatch::totalFlops() - jFlops;
00856   addFlops(flops + jFlops);
00857 
00858   SUNDANCE_MSG2(dfVerb(), 
00859     tab0 << "done QuadratureEvalMediator::fillFunctionCache()");
00860 }
00861 
00862 void QuadratureEvalMediator::computePhysQuadPts() const 
00863 {
00864   if (cacheIsValid()) 
00865   {
00866     Tabs tab0(0);
00867     SUNDANCE_MSG2(verb(), tab0 << "reusing phys quad points");
00868   }
00869   else
00870   {
00871     Tabs tab0(0);
00872     double jFlops = CellJacobianBatch::totalFlops();
00873     SUNDANCE_MSG2(verb(), tab0 << "computing phys quad points");
00874     physQuadPts_.resize(0);
00875     if (cellDim() != maxCellDim() && ElementIntegral::alwaysUseCofacets()
00876       && !isInternalBdry())
00877     {
00878       Tabs tab1;
00879       SUNDANCE_MSG2(verb(), tab1 << "using cofacets");
00880       SUNDANCE_MSG3(verb(), tab1 << "num cofacets = " << cofacetCellLID()->size());
00881       Array<Point> tmp;
00882       Array<int> cell(1);
00883       for (int c=0; c<cellLID()->size(); c++)
00884       {
00885         cell[0] = (*cofacetCellLID())[c];
00886         int facetIndex = facetIndices()[c];
00887         const Array<Point>& refFacetPts 
00888           =  (*(quadPtsReferredToMaxCell_.get(cellType())))[facetIndex];
00889         tmp.resize(refFacetPts.size());
00890         mesh().pushForward(maxCellDim(), cell, refFacetPts, tmp);
00891         for (int q=0; q<tmp.size(); q++)
00892         {
00893           physQuadPts_.append(tmp[q]);
00894         }
00895       }
00896     }
00897     else
00898     {
00899       const Array<Point>& refPts 
00900         = *(quadPtsForReferenceCell_.get(cellType()));
00901       mesh().pushForward(cellDim(), *cellLID(), 
00902         refPts, physQuadPts_); 
00903     }
00904     addFlops(CellJacobianBatch::totalFlops() - jFlops);
00905     cacheIsValid() = true;
00906     SUNDANCE_OUT(this->verb() > 2, 
00907       "phys quad: " << physQuadPts_);
00908   }
00909 }
00910 
00911 
00912 void QuadratureEvalMediator::print(std::ostream& os) const 
00913 {
00914   if (cacheIsValid())
00915   {
00916     Tabs tab0;
00917     os << tab0 << "Physical quadrature points" << std::endl;
00918     int i=0;
00919     for (int c=0; c<cellLID()->size(); c++)
00920     {
00921       Tabs tab1;
00922       os << tab1 << "cell " << c << std::endl;
00923       for (int q=0; q<physQuadPts_.size()/cellLID()->size(); q++, i++)
00924       {
00925         Tabs tab2;
00926         os << tab2 << "q=" << q << " " << physQuadPts_[i] << std::endl;
00927       }
00928     }
00929   }
00930 }
00931 
00932 

Site Contact