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 "SundanceCurveEvalMediator.hpp"
00043 #include "SundanceCoordExpr.hpp"
00044 #include "SundanceTempStack.hpp"
00045 #include "SundanceCellDiameterExpr.hpp"
00046 #include "SundanceCurveNormExpr.hpp"
00047 #include "SundanceCellVectorExpr.hpp"
00048 #include "SundanceSpatialDerivSpecifier.hpp"
00049 #include "SundanceDiscreteFunction.hpp"
00050 #include "SundanceDiscreteFuncElement.hpp"
00051 #include "SundanceCellJacobianBatch.hpp"
00052 #include "SundanceOut.hpp"
00053 #include "PlayaTabs.hpp"
00054 #include "PlayaExceptions.hpp"
00055 #include "SundanceCurveIntegralCalc.hpp"
00056
00057 #include "Teuchos_BLAS.hpp"
00058
00059
00060 using namespace Sundance;
00061 using namespace Teuchos;
00062 using namespace Playa;
00063
00064 TEUCHOS_TIMER(coordEvalTimer, "Quad mediator: coord eval")
00065
00066 CurveEvalMediator::CurveEvalMediator(
00067 const Mesh& mesh, const ParametrizedCurve& paramcurve ,
00068 int cellDim, const QuadratureFamily& quad) : StdFwkEvalMediator(mesh, cellDim) ,
00069 numEvaluationCases_(-1),
00070 quad_(quad),
00071 numQuadPtsForMaxCell_(0) ,
00072 paramcurve_(paramcurve) {
00073
00074
00075 maxCellType_ = mesh.cellType(mesh.spatialDim() );
00076 curveCellType_ = mesh.cellType(mesh.spatialDim()-1);
00077
00078
00079 Array<Point> quadPoints;
00080 Array<double> quadWeights;
00081 quad.getPoints(curveCellType_ , quadPoints , quadWeights );
00082 numQuadPtsForMaxCell_ = quadWeights.size();
00083 }
00084
00085 Time& CurveEvalMediator::coordEvaluationTimer()
00086 {
00087 return coordEvalTimer();
00088 }
00089
00090 void CurveEvalMediator::setCellType(const CellType& cellType,
00091 const CellType& maxCellType,
00092 bool isInternBdry)
00093 {
00094 StdFwkEvalMediator::setCellType(cellType, maxCellType, isInternBdry);
00095
00096 maxCellType_ = maxCellType;
00097
00098 Tabs tab;
00099
00100 SUNDANCE_MSG2(verb(), tab << "CurveEvalMediator::setCellType: cellType="
00101 << cellType << " cellDim=" << cellDim() << " maxCellType=" << maxCellType);
00102 SUNDANCE_MSG2(verb(), tab << "integration spec: ="
00103 << integrationCellSpec());
00104 SUNDANCE_MSG2(verb(), tab << "forbid cofacet integrations ="
00105 << forbidCofacetIntegrations());
00106 if (isInternalBdry())
00107 {
00108 SUNDANCE_MSG2(verb(), tab << "working on internal boundary");
00109 }
00110
00111
00112 if (cellType != maxCellType)
00113 {
00114
00115 }
00116 else
00117 {
00118 Tabs tab1;
00119 SUNDANCE_MSG2(verb(), tab1 << "no need for facet cases; work on original cell");
00120 numEvaluationCases_ = 1;
00121 }
00122
00123 SUNDANCE_MSG2(verb(), tab << "num eval cases =" << numEvaluationCases_);
00124
00125 }
00126
00127 int CurveEvalMediator::numQuadPts(const CellType& ct) const
00128 {
00129 if (ct != maxCellType_){
00130
00131 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "CurveEvalMediator::numQuadPts , cellType must be maxCellType_:" << ct);
00132 }
00133 return numQuadPtsForMaxCell_;
00134 }
00135
00136 void CurveEvalMediator::evalCellDiameterExpr(const CellDiameterExpr* expr,
00137 RCP<EvalVector>& vec) const
00138 {
00139
00140 Tabs tabs;
00141 SUNDANCE_MSG2(verb(),tabs
00142 << "CurveEvalMediator evaluating cell diameter expr "
00143 << expr->toString());
00144
00145 int nQuad = numQuadPts(cellType());
00146 int nCells = cellLID()->size();
00147
00148 SUNDANCE_MSG3(verb(),tabs << "number of quad pts=" << nQuad);
00149 Array<double> diameters;
00150 mesh().getCellDiameters(cellDim(), *cellLID(), diameters);
00151
00152 vec->resize(nQuad*nCells);
00153 double * const xx = vec->start();
00154 int k=0;
00155 for (int c=0; c<nCells; c++)
00156 {
00157 double h = diameters[c];
00158 for (int q=0; q<nQuad; q++, k++)
00159 {
00160 xx[k] = h;
00161 }
00162 }
00163 }
00164
00165 void CurveEvalMediator::evalCellVectorExpr(const CellVectorExpr* expr,
00166 RCP<EvalVector>& vec) const
00167 {
00168 Tabs tabs;
00169 SUNDANCE_MSG2(verb(),tabs
00170 << "CurveEvalMediator evaluating cell vector expr "
00171 << expr->toString());
00172
00173 int nQuad = numQuadPts(cellType());
00174 int nCells = cellLID()->size();
00175
00176 vec->resize(nQuad*nCells);
00177
00178 SUNDANCE_MSG3(verb(),tabs << "number of quad pts=" << nQuad);
00179 int dir = expr->componentIndex();
00180
00181 Array<Point> vectors;
00182 if (expr->isNormal())
00183 {
00184 mesh().outwardNormals(*cellLID(), vectors);
00185 }
00186 else
00187 {
00188 TEUCHOS_TEST_FOR_EXCEPTION(cellDim() != 1, std::runtime_error,
00189 "unable to compute tangent vectors for cell dim = " << cellDim());
00190 mesh().tangentsToEdges(*cellLID(), vectors);
00191 }
00192
00193 double * const xx = vec->start();
00194 int k=0;
00195 for (int c=0; c<nCells; c++)
00196 {
00197 double n = vectors[c][dir];
00198 for (int q=0; q<nQuad; q++, k++)
00199 {
00200 xx[k] = n;
00201 }
00202 }
00203 }
00204
00205 void CurveEvalMediator::evalCoordExpr(const CoordExpr* expr,
00206 RCP<EvalVector>& vec) const
00207 {
00208 Tabs tabs;
00209 SUNDANCE_MSG2(verb(),tabs
00210 << "CurveEvalMediator evaluating coord expr "
00211 << expr->toString());
00212
00213 TimeMonitor timer(coordEvalTimer());
00214
00215 int nQuad = numQuadPts(cellType());
00216 int nCells = cellLID()->size();
00217 int d = expr->dir();
00218
00219 SUNDANCE_MSG3(verb(),tabs << "number of quad pts=" << nQuad << ", nCells:" << nCells);
00220
00221
00222 vec->resize(nQuad*nCells);
00223 double * const xx = vec->start();
00224 int k=0;
00225 int maxDim = mesh().spatialDim() ;
00226 Array<int> cellLIDs_tmp(1);
00227 Array<Point> phyPoints;
00228 Array<Point> refPoints;
00229 Array<Point> refDevs;
00230 Array<Point> refNormal;
00231
00232
00233 for (int c=0; c<nCells; c++)
00234 {
00235 int maxCellLID = (*cellLID())[c];
00236 cellLIDs_tmp[0] = (*cellLID())[c];
00237
00238 if ( mesh().hasCurvePoints( maxCellLID , paramcurve_.myID() ))
00239 {
00240 mesh().getCurvePoints( maxCellLID , paramcurve_.myID() , refPoints , refDevs , refNormal );
00241 }
00242 else
00243 {
00244
00245 CurveIntegralCalc::getCurveQuadPoints( maxCellType_ , maxCellLID , mesh() , paramcurve_ , quad_ ,
00246 refPoints, refDevs , refNormal);
00247
00248 mesh().setCurvePoints( maxCellLID , paramcurve_.myID() , refPoints, refDevs , refNormal );
00249 }
00250
00251
00252 mesh().pushForward( maxDim , cellLIDs_tmp , refPoints , phyPoints );
00253
00254
00255 for (int q=0; q<nQuad; q++, k++)
00256 {
00257 xx[k] = phyPoints[q][d];
00258 }
00259 }
00260 }
00261
00262
00263 void CurveEvalMediator::evalCurveNormExpr(const CurveNormExpr* expr,
00264 RCP<EvalVector>& vec) const {
00265
00266 Tabs tabs;
00267 SUNDANCE_MSG2(verb(),tabs
00268 << "CurveEvalMediator evaluating curve normal expr "
00269 << expr->toString());
00270
00271 int nQuad = numQuadPts(cellType());
00272 int nCells = cellLID()->size();
00273 int d = expr->dir();
00274
00275 SUNDANCE_MSG3(verb(),tabs << "number of quad pts=" << nQuad);
00276
00277
00278 vec->resize(nCells*nQuad);
00279
00280 double * const xx = vec->start();
00281 int k=0;
00282 Array<int> cellLIDs_tmp(1);
00283 Array<Point> phyPoints;
00284 Array<Point> refPoints;
00285 Array<Point> refDevs;
00286 Array<Point> refNormal;
00287
00288
00289 for (int c=0; c<nCells; c++)
00290 {
00291 int maxCellLID = (*cellLID())[c];
00292 cellLIDs_tmp[0] = (*cellLID())[c];
00293
00294 if ( mesh().hasCurvePoints( maxCellLID , paramcurve_.myID() ))
00295 {
00296 mesh().getCurvePoints( maxCellLID , paramcurve_.myID() , refPoints , refDevs , refNormal );
00297 }
00298 else
00299 {
00300
00301 CurveIntegralCalc::getCurveQuadPoints( maxCellType_ , maxCellLID , mesh() , paramcurve_ , quad_ ,
00302 refPoints, refDevs , refNormal);
00303
00304 mesh().setCurvePoints( maxCellLID , paramcurve_.myID() , refPoints, refDevs , refNormal );
00305 }
00306
00307
00308 for (int q=0; q<nQuad; q++, k++)
00309 {
00310 xx[k] = refNormal[q][d];
00311 }
00312 }
00313 }
00314
00315
00316 void CurveEvalMediator
00317 ::evalDiscreteFuncElement(const DiscreteFuncElement* expr,
00318 const Array<MultiIndex>& multiIndices,
00319 Array<RCP<EvalVector> >& vec) const
00320 {
00321
00322 int verbo = dfVerb();
00323 Tabs tab1;
00324 SUNDANCE_MSG2(verbo , tab1
00325 << "CurveEvalMediator evaluating Discrete Function expr " << expr->toString());
00326
00327 const DiscreteFunctionData* f = DiscreteFunctionData::getData(expr);
00328
00329 TEUCHOS_TEST_FOR_EXCEPTION(f==0, std::logic_error,
00330 "QuadratureEvalMediator::evalDiscreteFuncElement() called "
00331 "with expr that is not a discrete function");
00332
00333 SUNDANCE_MSG2(verbo , tab1 << "After casting DiscreteFunctionData" << expr->toString());
00334
00335 RCP<Array<Array<double> > > localValues;
00336 Array<int> cellLIDs_tmp(1);
00337 Array<Point> phyPoints;
00338 Array<Point> refPoints;
00339 Array<Point> refDevs;
00340 Array<Point> refNormal;
00341 int nCells = cellLID()->size();
00342
00343 RCP<const MapStructure> mapStruct;
00344 int myIndex = expr->myIndex();
00345 int nQuad = numQuadPtsForMaxCell_;
00346 Array<int> k(multiIndices.size(),0);
00347
00348 Teuchos::BLAS<int,double> blas;
00349
00350 SUNDANCE_MSG2(verbo , tab1 << "After declaring BLAS: " << expr->toString());
00351
00352
00353 for (int i=0; i<multiIndices.size(); i++)
00354 {
00355 vec[i]->resize(nCells*nQuad);
00356 }
00357
00358
00359 for (int c=0; c<nCells; c++)
00360 {
00361 int maxCellLID = (*cellLID())[c];
00362 localValues = rcp(new Array<Array<double> >());
00363 SUNDANCE_MSG2(verbo , tab1 << "Cell:" << c << " of " << nCells << " , maxCellLID:" << maxCellLID );
00364
00365 cellLIDs_tmp[0] = (*cellLID())[c];
00366 SUNDANCE_MSG2(verbo , tab1 << " Before calling f->getLocalValues:" << cellLIDs_tmp.size() <<
00367 " f==0 : " << (f==0) << " tmp:" << f->mesh().spatialDim());
00368
00369 mapStruct = f->getLocalValues(maxCellDim(), cellLIDs_tmp , *localValues);
00370 SUNDANCE_MSG2(verbo , tab1 << " After getting mapStruct:" << maxCellLID );
00371 SUNDANCE_MSG2(verbo , tab1 << " mapStruct->numBasisChunks():" << mapStruct->numBasisChunks() );
00372 int chunk = mapStruct->chunkForFuncID(myIndex);
00373 int funcIndex = mapStruct->indexForFuncID(myIndex);
00374 int nFuncs = mapStruct->numFuncs(chunk);
00375 SUNDANCE_MSG2(verbo , tab1 << " chunk:" << chunk );
00376 SUNDANCE_MSG2(verbo , tab1 << " funcIndex:" << funcIndex );
00377 SUNDANCE_MSG2(verbo , tab1 << " nFuncs:" << nFuncs );
00378
00379
00380 BasisFamily basis = rcp_dynamic_cast<BasisFamilyBase>(mapStruct->basis(chunk));
00381 int nNodesTotal = basis.nReferenceDOFsWithFacets(maxCellType(), maxCellType());
00382
00383
00384 if ( mesh().hasCurvePoints( maxCellLID , paramcurve_.myID() ))
00385 {
00386 mesh().getCurvePoints( maxCellLID , paramcurve_.myID() , refPoints , refDevs , refNormal );
00387 }
00388 else
00389 {
00390
00391 CurveIntegralCalc::getCurveQuadPoints( maxCellType_ , maxCellLID , mesh() , paramcurve_ , quad_ ,
00392 refPoints, refDevs , refNormal);
00393
00394 mesh().setCurvePoints( maxCellLID , paramcurve_.myID() , refPoints, refDevs , refNormal );
00395 }
00396
00397
00398 SUNDANCE_MSG2(verbo , tab1 << " multiIndices.size()" << multiIndices.size() );
00399 for (int i=0; i<multiIndices.size(); i++)
00400 {
00401
00402 int nDerivResults = 1;
00403 if ( multiIndices[i].order() == 1 ) nDerivResults = maxCellDim();
00404
00405 int pDir = 0;
00406 int derivNum = 1;
00407 MultiIndex mi;
00408 SUNDANCE_MSG2(verbo , tab1 << " before asking anything i = " << i);
00409 SUNDANCE_MSG2(verbo , tab1 << " multiindex order : " << multiIndices[i].order());
00410 SUNDANCE_MSG2(verbo , tab1 << " multiindex : " << multiIndices[i] );
00411 if (multiIndices[i].order() > 0){
00412 pDir = multiIndices[i].firstOrderDirection();
00413 mi[pDir] = 1;
00414 derivNum = mesh().spatialDim();
00415 }
00416
00417 Array<Array<double> > result(nQuad*derivNum);
00418 Array<Array<Array<double> > > tmp;
00419
00420 int offs = nNodesTotal;
00421
00422
00423 for (int deriv = 0 ; deriv < derivNum ; deriv++)
00424 {
00425
00426 if (multiIndices[i].order() > 0){
00427
00428 MultiIndex mi_tmp;
00429 mi_tmp[deriv] = 1;
00430 SpatialDerivSpecifier deriv(mi_tmp);
00431 SUNDANCE_MSG2(verbo , tab1 << "computing derivatives : " << deriv << " on reference cell ");
00432 basis.refEval( maxCellType_ , refPoints , deriv, tmp , verbo );
00433 } else {
00434 SpatialDerivSpecifier deriv(mi);
00435
00436 SUNDANCE_MSG2(verbo , tab1 << "computing values reference cell ");
00437 basis.refEval( maxCellType_ , refPoints , deriv, tmp , verbo );
00438 }
00439
00440 SUNDANCE_MSG2(verbo , tab1 << "resize result vector , offs:" << offs);
00441 for (int q=0; q<nQuad; q++){
00442 result[nQuad*deriv + q].resize(offs);
00443 }
00444
00445
00446 SUNDANCE_MSG2(verbo , tab1 << "copy results ");
00447 int offs1 = 0;
00448 for (int q=0; q<nQuad; q++)
00449 {
00450 offs1 = 0;
00451 for (int d=0; d<basis.dim(); d++)
00452 {
00453 int nNodes = tmp[d][q].size();
00454 for (int n=0; n<nNodes; n++ , offs1++ )
00455 {
00456 result[nQuad*deriv + q][offs1] = tmp[d][q][n];
00457 }
00458 }
00459 }
00460 }
00461
00462
00463 SUNDANCE_MSG2(verbo , tab1 << "summing up values , funcIndex:" << funcIndex << " offs:" << offs);
00464 for (int deriv = 0 ; deriv < derivNum ; deriv++)
00465 {
00466 for (int q=0; q<nQuad; q++)
00467 {
00468 double sum = 0.0;
00469
00470 for (int n = 0 ; n < offs ; n++){
00471 sum = sum + result[nQuad*deriv + q][n] * (*localValues)[chunk][funcIndex*offs + n];
00472 }
00473
00474 result[nQuad*deriv + q][0] = sum;
00475 }
00476 }
00477
00478
00479 const CellJacobianBatch& J = JTrans();
00480 if (mi.order()==1)
00481 {
00482 Tabs tab1;
00483 Tabs tab2;
00484 SUNDANCE_MSG2(verbo, tab2 << "Jacobian batch nCells=" << J.numCells());
00485 SUNDANCE_MSG2(verbo, tab2 << "Jacobian batch cell dim=" << J.cellDim());
00486 SUNDANCE_MSG2(verbo, tab2 << "Jacobian batch spatial dim=" << J.spatialDim());
00487
00488 Array<double> invJ;
00489 J.getInvJ(c, invJ);
00490 for (int q=0; q<nQuad; q++)
00491 {
00492 double sum = 0.0;
00493 for (int deriv = 0 ; deriv < derivNum ; deriv++)
00494 {
00495
00496 sum = sum + result[nQuad*deriv + q][0] * invJ[derivNum*pDir + deriv];
00497 }
00498
00499 result[q][0] = sum;
00500 }
00501 }
00502
00503
00504
00505 double* vecPtr = vec[i]->start();
00506 for (int q=0; q<nQuad; q++, k[i]++)
00507 {
00508 vecPtr[k[i]] = result[q][0];
00509 }
00510 SUNDANCE_MSG2(verbo , tab1 << " END copy results back ");
00511 }
00512 SUNDANCE_MSG2(verbo , tab1 << " END loop over multiindex ");
00513 }
00514 SUNDANCE_MSG2(verbo , tab1 << " END loop over cells ");
00515 }
00516
00517
00518 void CurveEvalMediator::print(std::ostream& os) const
00519 {
00520
00521 }
00522
00523