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