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 #include "SundanceCurveQuadratureIntegral.hpp"
00032
00033 #include "SundanceCurveIntegralCalc.hpp"
00034 #include "SundanceGaussianQuadrature.hpp"
00035 #include "SundanceSpatialDerivSpecifier.hpp"
00036 #include "SundanceOut.hpp"
00037 #include "SundanceTabs.hpp"
00038 #include "Teuchos_TimeMonitor.hpp"
00039
00040 using namespace Sundance;
00041 using namespace Teuchos;
00042
00043 using std::endl;
00044 using std::setw;
00045 using std::setprecision;
00046
00047 extern "C"
00048 {
00049 int dgemm_(const char* transA, const char* transB,
00050 const int* M, const int *N, const int* K,
00051 const double* alpha,
00052 const double* A, const int* ldA,
00053 const double* B, const int* ldB,
00054 const double* beta,
00055 double* C, const int* ldC);
00056 }
00057
00058 static Time& maxCellQuadratureTimer()
00059 {
00060 static RCP<Time> rtn
00061 = TimeMonitor::getNewTimer("max cell quadrature");
00062 return *rtn;
00063 }
00064
00065
00066 CurveQuadratureIntegral::CurveQuadratureIntegral(
00067 const CellType& cellType,
00068 const bool isConstantIntegral,
00069 const QuadratureFamily& quad,
00070 const ParametrizedCurve& globalCurve,
00071 const Mesh& mesh,
00072 int verb)
00073 : ElementIntegral(dimension(cellType), cellType, dimension(cellType),
00074 cellType, true, globalCurve, mesh, verb),
00075 quad_(quad),
00076 quadPts_(),
00077 quadWeights_(),
00078 W_(),
00079 useSumFirstMethod_(true),
00080 useConstCoeff_(isConstantIntegral)
00081 {
00082 Tabs tab0(0);
00083
00084 TEST_FOR_EXCEPTION( cellType != mesh.cellType(mesh.spatialDim()) , RuntimeError, " CurveQuadratureIntegral::CurveQuadratureIntegral , only on MAXcell!!");
00085
00086 SUNDANCE_MSG1(setupVerb(), tab0 << "CurveQuadratureIntegral ctor for 0-form");
00087 if (setupVerb()) describe(Out::os());
00088
00089 SUNDANCE_MSG1(setupVerb(), tab0 << "quadrature family=" << quad);
00090
00091
00092 quad.getPoints(mesh.cellType(mesh.spatialDim()-1), quadPts_, quadWeights_);
00093 int nQuad = quadPts_.size();
00094
00095 W_.resize(nQuad);
00096 quadCurveDerivs_.resize(nQuad);
00097 quadCurveNormals_.resize(nQuad);
00098
00099 for (int q=0; q<nQuad; q++)
00100 {
00101 W_[q] = quadWeights_[q];
00102 }
00103 }
00104
00105
00106 CurveQuadratureIntegral::CurveQuadratureIntegral(
00107 const CellType& cellType,
00108 const bool isConstantIntegral,
00109 const BasisFamily& testBasis,
00110 int alpha,
00111 int testDerivOrder,
00112 const QuadratureFamily& quad,
00113 const ParametrizedCurve& globalCurve,
00114 const Mesh& mesh,
00115 int verb)
00116 : ElementIntegral(dimension(cellType), cellType, dimension(cellType),
00117 cellType, testBasis,
00118 alpha, testDerivOrder, true, globalCurve, mesh, verb),
00119 quad_(quad),
00120 quadPts_(),
00121 quadWeights_(),
00122 W_(),
00123 useSumFirstMethod_(true),
00124 useConstCoeff_(isConstantIntegral)
00125 {
00126 Tabs tab0(0);
00127
00128 TEST_FOR_EXCEPTION( cellType != mesh.cellType(mesh.spatialDim()) , RuntimeError, " CurveQuadratureIntegral::CurveQuadratureIntegral , only on MAXcell!!");
00129
00130 SUNDANCE_MSG1(setupVerb(), tab0 << "CurveQuadratureIntegral ctor for 1-form");
00131 if (setupVerb()) describe(Out::os());
00132 assertLinearForm();
00133
00134
00135 SUNDANCE_MSG1(setupVerb(), tab0 << "quadrature family=" << quad);
00136
00137 quad.getPoints(mesh.cellType(mesh.spatialDim()-1), quadPts_, quadWeights_);
00138 int nQuad = quadPts_.size();
00139
00140 W_.resize(nQuad * nRefDerivTest() * nNodesTest());
00141 quadCurveDerivs_.resize(nQuad);
00142 quadCurveNormals_.resize(nQuad);
00143
00144 SUNDANCE_MSG1(setupVerb(), tab0 << "num nodes for test function " << nNodesTest());
00145 }
00146
00147
00148
00149
00150 CurveQuadratureIntegral::CurveQuadratureIntegral(
00151 const CellType& cellType,
00152 const bool isConstantIntegral,
00153 const BasisFamily& testBasis,
00154 int alpha,
00155 int testDerivOrder,
00156 const BasisFamily& unkBasis,
00157 int beta,
00158 int unkDerivOrder,
00159 const QuadratureFamily& quad,
00160 const ParametrizedCurve& globalCurve,
00161 const Mesh& mesh,
00162 int verb)
00163 : ElementIntegral(dimension(cellType), cellType, dimension(cellType),
00164 cellType, testBasis,
00165 alpha, testDerivOrder, unkBasis, beta, unkDerivOrder, true,
00166 globalCurve, mesh,
00167 verb),
00168 quad_(quad),
00169 quadPts_(),
00170 quadWeights_(),
00171 W_(),
00172 useSumFirstMethod_(true),
00173 useConstCoeff_(isConstantIntegral)
00174 {
00175 Tabs tab0(0);
00176
00177 TEST_FOR_EXCEPTION( cellType != mesh.cellType(mesh.spatialDim()) , RuntimeError, " CurveQuadratureIntegral::CurveQuadratureIntegral , only on MAXcell!!");
00178
00179 SUNDANCE_MSG1(setupVerb(), tab0 << "CurveQuadratureIntegral ctor for 2-form");
00180 if (setupVerb()) describe(Out::os());
00181 assertBilinearForm();
00182
00183
00184 quad.getPoints(mesh.cellType(mesh.spatialDim()-1), quadPts_, quadWeights_);
00185 int nQuad = quadPts_.size();
00186
00187 W_.resize(nQuad * nRefDerivTest() * nNodesTest()
00188 * nRefDerivUnk() * nNodesUnk());
00189
00190 quadCurveDerivs_.resize(nQuad);
00191 quadCurveNormals_.resize(nQuad);
00192 }
00193
00194
00195
00196 void CurveQuadratureIntegral::updateRefCellInformation(int maxCellLID , const ParametrizedCurve& curve) const{
00197
00198
00199 Tabs tabs;
00200 SUNDANCE_MSG3(integrationVerb(), tabs << "CurveQuadratureIntegral::updateRefCellInformation, "
00201 "Update Cell LID: " << maxCellLID);
00202
00203 if ( mesh().hasCurvePoints( maxCellLID , curve.myID() ))
00204 {
00205 SUNDANCE_MSG3(integrationVerb(), tabs << "CurveQuadratureIntegral::updateRefCellInformation, has Points, just returning them");
00206 mesh().getCurvePoints( maxCellLID , curve.myID() , quadPts_ , quadCurveDerivs_ , quadCurveNormals_ );
00207 }
00208 else
00209 {
00210
00211 SUNDANCE_MSG3(integrationVerb(), tabs << "CurveQuadratureIntegral::updateRefCellInformation, computing the intersection points");
00212 CurveIntegralCalc::getCurveQuadPoints( mesh().cellType(mesh().spatialDim()) , maxCellLID , mesh() , curve ,
00213 quad_ , quadPts_ , quadCurveDerivs_ , quadCurveNormals_ );
00214
00215 SUNDANCE_MSG3(integrationVerb(), tabs << "CurveQuadratureIntegral::updateRefCellInformation, storing the intersection points in the mesh");
00216 mesh().setCurvePoints( maxCellLID , curve.myID() , quadPts_ , quadCurveDerivs_ , quadCurveNormals_ );
00217 }
00218 }
00219
00220
00221 void CurveQuadratureIntegral::updateRefCellIntegralOneForm(int maxCellLID , int cellInBatch) const {
00222
00223
00224 Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());
00225
00226 Tabs tabs;
00227 SUNDANCE_MSG2(integrationVerb(), tabs << "CurveQuadratureIntegral::updateRefCellIntegralTwoForm, "
00228 "Update Cell LID: " << maxCellLID);
00229
00230
00231 updateRefCellInformation( maxCellLID , globalCurve() );
00232
00233 for (int r=0; r<nRefDerivTest(); r++)
00234 {
00235 Tabs tab3;
00236 SUNDANCE_MSG1(setupVerb(), tab3
00237 << "evaluating basis functions for ref deriv direction " << r);
00238 MultiIndex mi;
00239 testBasisVals[r].resize(testBasis().dim());
00240 if (testDerivOrder()==1) mi[r] = 1;
00241 SpatialDerivSpecifier deriv(mi);
00242 testBasis().refEval(evalCellType(), quadPts_, deriv,
00243 testBasisVals[r], setupVerb());
00244 }
00245
00246 int vecComp = 0;
00247 int nQuad = quadPts_.size();
00248 for (int q=0; q<nQuad; q++)
00249 {
00250 for (int t=0; t<nRefDerivTest(); t++)
00251 {
00252 for (int nt=0; nt<nNodesTest(); nt++)
00253 {
00254 wValue(q, t, nt)
00255 = chop(quadWeights_[q] * testBasisVals[t][vecComp][q][nt]) ;
00256 }
00257 }
00258 }
00259
00260
00261 for (int i=0; i<W_.size(); i++) W_[i] = chop(W_[i]);
00262 }
00263
00264
00265 void CurveQuadratureIntegral::updateRefCellIntegralTwoForm(int maxCellLID , int cellInBatch) const {
00266
00267
00268 Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());
00269 Array<Array<Array<Array<double> > > > unkBasisVals(nRefDerivUnk());
00270
00271 Tabs tabs;
00272 SUNDANCE_MSG2(integrationVerb(), tabs << "CurveQuadratureIntegral::updateRefCellIntegralTwoForm, "
00273 "Update Cell LID: " << maxCellLID);
00274
00275
00276 updateRefCellInformation( maxCellLID , globalCurve() );
00277
00278 for (int r=0; r<nRefDerivTest(); r++)
00279 {
00280 testBasisVals[r].resize(testBasis().dim());
00281 MultiIndex mi;
00282 if (testDerivOrder()==1) mi[r] = 1;
00283 SpatialDerivSpecifier deriv(mi);
00284 testBasis().refEval(evalCellType(), quadPts_, deriv,
00285 testBasisVals[r], setupVerb());
00286 }
00287
00288 for (int r=0; r<nRefDerivUnk(); r++)
00289 {
00290 unkBasisVals[r].resize(unkBasis().dim());
00291 MultiIndex mi;
00292 if (unkDerivOrder()==1) mi[r] = 1;
00293 SpatialDerivSpecifier deriv(mi);
00294 unkBasis().refEval(evalCellType(),
00295 quadPts_, deriv, unkBasisVals[r], setupVerb());
00296 }
00297
00298
00299 int vecComp = 0;
00300 int nQuad = quadPts_.size();
00301
00302 for (int q=0; q<nQuad; q++)
00303 {
00304 for (int t=0; t<nRefDerivTest(); t++)
00305 {
00306 for (int nt=0; nt<nNodesTest(); nt++)
00307 {
00308 for (int u=0; u<nRefDerivUnk(); u++)
00309 {
00310 for (int nu=0; nu<nNodesUnk(); nu++)
00311 {
00312 wValue(q, t, nt, u, nu)
00313 = chop(quadWeights_[q] * testBasisVals[t][vecComp][q][nt]
00314 * unkBasisVals[u][vecComp][q][nu]);
00315 }
00316 }
00317 }
00318 }
00319 }
00320
00321
00322 for (int i=0; i<W_.size(); i++) W_[i] = chop(W_[i]);
00323 }
00324
00325 void CurveQuadratureIntegral
00326 ::transformZeroForm(const CellJacobianBatch& JTrans,
00327 const CellJacobianBatch& JVol,
00328 const Array<int>& isLocalFlag,
00329 const Array<int>& facetIndex,
00330 const RCP<Array<int> >& cellLIDs,
00331 const double constCoeff,
00332 const double* const coeff,
00333 RCP<Array<double> >& A) const
00334 {
00335 TimeMonitor timer(maxCellQuadratureTimer());
00336 Tabs tabs;
00337 SUNDANCE_MSG1(integrationVerb(), tabs << "doing zero form by quadrature");
00338
00339 TEST_FOR_EXCEPTION(order() != 0, InternalError,
00340 "CurveQuadratureIntegral::transformZeroForm() called "
00341 "for form of order " << order());
00342
00343 TEST_FOR_EXCEPTION( (int) isLocalFlag.size() != 0
00344 && (int) isLocalFlag.size() != JVol.numCells(),
00345 RuntimeError,
00346 "mismatch between isLocalFlag.size()=" << isLocalFlag.size()
00347 << " and JVol.numCells()=" << JVol.numCells());
00348
00349 bool checkLocalFlag = (int) isLocalFlag.size() != 0;
00350
00351 int nQuad = quadWeights_.size();
00352
00353
00354 double* coeffPtr = (double*) coeff;
00355 Array<double> constCoeff_arr( JVol.numCells() * nQuad);
00356 if (useConstCoeff_){
00357 for (int jk=0; jk < JVol.numCells() * nQuad ; jk++ ) { constCoeff_arr[jk] = constCoeff; };
00358 coeffPtr = &(constCoeff_arr[0]);
00359 }
00360
00361 double& a = (*A)[0];
00362 SUNDANCE_MSG5(integrationVerb(), tabs << "input A=");
00363 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00364 const Array<double>& w = W_;
00365
00366
00367 for (int c=0; c<JVol.numCells(); c++)
00368 {
00369
00370 if (checkLocalFlag && !isLocalFlag[c])
00371 {
00372 coeffPtr += nQuad;
00373 continue;
00374 }
00375
00376
00377
00378 for (int q=0; q<nQuad; q++, coeffPtr++)
00379 {
00380
00381
00382 a += w[q]*(*coeffPtr);
00383 }
00384 }
00385
00386 SUNDANCE_MSG5(integrationVerb(), tabs << "output A = ");
00387 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00388
00389 SUNDANCE_MSG1(integrationVerb(), tabs << "done zero form");
00390 }
00391
00392
00393 void CurveQuadratureIntegral::transformOneForm(const CellJacobianBatch& JTrans,
00394 const CellJacobianBatch& JVol,
00395 const Array<int>& facetIndex,
00396 const RCP<Array<int> >& cellLIDs,
00397 const double constCoeff,
00398 const double* const coeff,
00399 RCP<Array<double> >& A) const
00400 {
00401 TimeMonitor timer(maxCellQuadratureTimer());
00402 Tabs tabs;
00403 TEST_FOR_EXCEPTION(order() != 1, InternalError,
00404 "CurveQuadratureIntegral::transformOneForm() called for form "
00405 "of order " << order());
00406 SUNDANCE_MSG2(integrationVerb(), tabs << "doing one form by quadrature");
00407 int flops = 0;
00408
00409 int nQuad = quadWeights_.size();
00410
00411
00412
00413
00414 if (testDerivOrder() == 0)
00415 {
00416 double* aPtr = &((*A)[0]);
00417 SUNDANCE_MSG5(integrationVerb(), tabs << "input A = ");
00418 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00419
00420 int offset = 0 ;
00421 const Array<double>& w = W_;
00422 double* coeffPtr = (double*) coeff;
00423
00424 Array<double> constCoeff_arr( JVol.numCells() * nQuad);
00425 if (useConstCoeff_){
00426 for (int jk=0; jk < JVol.numCells() * nQuad ; jk++ ) { constCoeff_arr[jk] = constCoeff; };
00427 coeffPtr = &(constCoeff_arr[0]);
00428 }
00429
00430 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00431 {
00432 Tabs tab2;
00433 SUNDANCE_MSG4(integrationVerb(), tab2 << "c=" << c << " detJ=" << JVol.detJ()[c] );
00434
00435
00436 updateRefCellIntegralOneForm( (*cellLIDs)[c] , c );
00437
00438 for (int q=0; q<nQuad; q++, coeffPtr++)
00439 {
00440 Tabs tab3;
00441
00442 double curveDerivNorm = sqrt(quadCurveDerivs_[q]*quadCurveDerivs_[q]);
00443 double f = (*coeffPtr) * curveDerivNorm;
00444 SUNDANCE_MSG4(integrationVerb(), tab3 << "q=" << q << " coeff=" <<
00445 *coeffPtr << " coeff*detJ=" << f);
00446 for (int n=0; n<nNodes(); n++)
00447 {
00448 Tabs tab4;
00449 SUNDANCE_MSG4(integrationVerb(), tab4 << "n=" << n << " w=" <<
00450 w[n + nNodes()*q]);
00451 aPtr[offset+n] += f*w[n + nNodes()*q];
00452 }
00453 }
00454
00455 if (integrationVerb() >= 4)
00456 {
00457 Out::os() << tab2 << "integration results on cell:" << std::endl;
00458 Out::os() << tab2 << setw(10) << "n" << setw(12) << "I_n" << std::endl;
00459 for (int n=0; n<nNodes(); n++)
00460 {
00461 Out::os() << tab2 << setw(10) << n
00462 << setw(12) << setprecision(6) << aPtr[offset+n] << std::endl;
00463 }
00464 }
00465
00466 }
00467
00468 SUNDANCE_MSG5(integrationVerb(), tabs << "output A = ");
00469 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00470 }
00471 else
00472 {
00473
00474
00475 createOneFormTransformationMatrix(JTrans, JVol);
00476
00477 SUNDANCE_MSG4(transformVerb(),
00478 Tabs() << "transformation matrix=" << G(alpha()));
00479
00480 double* GPtr = &(G(alpha())[0]);
00481
00482 transformSummingFirst(JVol.numCells(), JVol ,facetIndex, cellLIDs, constCoeff , GPtr, coeff, A);
00483
00484 }
00485 addFlops(flops);
00486 }
00487
00488
00489 void CurveQuadratureIntegral::transformTwoForm(const CellJacobianBatch& JTrans,
00490 const CellJacobianBatch& JVol,
00491 const Array<int>& facetIndex,
00492 const RCP<Array<int> >& cellLIDs,
00493 const double constCoeff,
00494 const double* const coeff,
00495 RCP<Array<double> >& A) const
00496 {
00497 TimeMonitor timer(maxCellQuadratureTimer());
00498 Tabs tabs;
00499 TEST_FOR_EXCEPTION(order() != 2, InternalError,
00500 "CurveQuadratureIntegral::transformTwoForm() called for form "
00501 "of order " << order());
00502 SUNDANCE_MSG2(integrationVerb(), tabs << "doing one form by quadrature");
00503
00504 int nQuad = quadWeights_.size();
00505
00506
00507
00508
00509 if (testDerivOrder() == 0 && unkDerivOrder() == 0)
00510 {
00511 double* aPtr = &((*A)[0]);
00512 int offset = 0 ;
00513 const Array<double>& w = W_;
00514 double* coeffPtr = (double*) coeff;
00515
00516
00517 Array<double> constCoeff_arr( JVol.numCells() * nQuad);
00518 SUNDANCE_MSG2(integrationVerb(), tabs << "CurveQuadratureIntegral::transformTwoForm , use const variabl: " << useConstCoeff_ );
00519 if (useConstCoeff_){
00520 SUNDANCE_MSG2(integrationVerb(), tabs << "CurveQuadratureIntegral::transformTwoForm , "
00521 "Using constant coeff: " << constCoeff);
00522 for (int jk=0; jk < JVol.numCells() * nQuad ; jk++ ) { constCoeff_arr[jk] = constCoeff; };
00523 coeffPtr = &(constCoeff_arr[0]);
00524 }
00525
00526 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00527 {
00528
00529 updateRefCellIntegralTwoForm( (*cellLIDs)[c] , c );
00530
00531 for (int q=0; q<nQuad; q++, coeffPtr++)
00532 {
00533
00534 double curveDerivNorm = sqrt(quadCurveDerivs_[q]*quadCurveDerivs_[q]);
00535 double f = (*coeffPtr)*curveDerivNorm;
00536 for (int n=0; n<nNodes(); n++)
00537 {
00538 aPtr[offset+n] += f*w[n + nNodes()*q];
00539 }
00540 }
00541 }
00542
00543 }
00544 else
00545 {
00546 createTwoFormTransformationMatrix(JTrans, JVol);
00547 double* GPtr;
00548
00549 if (testDerivOrder() == 0)
00550 {
00551 GPtr = &(G(beta())[0]);
00552 SUNDANCE_MSG2(transformVerb(),
00553 Tabs() << "transformation matrix=" << G(beta()));
00554 }
00555 else if (unkDerivOrder() == 0)
00556 {
00557 GPtr = &(G(alpha())[0]);
00558 SUNDANCE_MSG2(transformVerb(),
00559 Tabs() << "transformation matrix=" << G(alpha()));
00560 }
00561 else
00562 {
00563 GPtr = &(G(alpha(), beta())[0]);
00564 SUNDANCE_MSG2(transformVerb(),
00565 Tabs() << "transformation matrix="
00566 << G(alpha(),beta()));
00567 }
00568
00569 transformSummingFirst(JTrans.numCells(), JVol , facetIndex, cellLIDs, constCoeff , GPtr, coeff, A);
00570
00571 }
00572 }
00573
00574 void CurveQuadratureIntegral
00575 ::transformSummingFirst(int nCells,
00576 const CellJacobianBatch& JVol,
00577 const Array<int>& facetIndex,
00578 const RCP<Array<int> >& cellLIDs,
00579 const double constCoeff,
00580 const double* const GPtr,
00581 const double* const coeff,
00582 RCP<Array<double> >& A) const
00583 {
00584 double* aPtr = &((*A)[0]);
00585 int nQuad = quadWeights_.size();
00586 double* coeffPtr = (double*) coeff;
00587
00588
00589 Array<double> constCoeff_arr( JVol.numCells() * nQuad);
00590 if (useConstCoeff_){
00591 for (int jk=0; jk < JVol.numCells() * nQuad ; jk++ ) { constCoeff_arr[jk] = constCoeff; };
00592 coeffPtr = &(constCoeff_arr[0]);
00593 }
00594
00595 int transSize = 0;
00596 if (order()==2)
00597 {
00598 transSize = nRefDerivTest() * nRefDerivUnk();
00599 }
00600 else
00601 {
00602 transSize = nRefDerivTest();
00603 }
00604
00605
00606 static Array<double> sumWorkspace;
00607
00608 int swSize = transSize * nNodes();
00609 sumWorkspace.resize(swSize);
00610 const Array<double>& w = W_;
00611
00612
00613
00614
00615
00616
00617
00618
00619
00620
00621 for (int c=0; c<nCells; c++)
00622 {
00623
00624 for (int i=0; i<swSize; i++) sumWorkspace[i]=0.0;
00625
00626
00627 if (order() == 2) {
00628 updateRefCellIntegralTwoForm( (*cellLIDs)[c] , c );
00629 }
00630 if (order() == 1) {
00631 updateRefCellIntegralOneForm( (*cellLIDs)[c] , c );
00632 }
00633
00634 for (int q=0; q<nQuad; q++, coeffPtr++)
00635 {
00636
00637 double curveDerivNorm = sqrt(quadCurveDerivs_[q]*quadCurveDerivs_[q]);
00638
00639 double f = (*coeffPtr);
00640 for (int n=0; n<swSize; n++)
00641 {
00642 sumWorkspace[n] += f*w[n + q*swSize] * curveDerivNorm;
00643 }
00644 }
00645
00646
00647 const double * const gCell = &(GPtr[transSize*c]);
00648 double* aCell = aPtr + nNodes()*c;
00649
00650 double detJInv = 1/fabs(JVol.detJ()[c]);
00651 for (int i=0; i<nNodes(); i++)
00652 {
00653 for (int j=0; j<transSize; j++)
00654 {
00655 aCell[i] += sumWorkspace[nNodes()*j + i] * gCell[j] * detJInv;
00656 }
00657 }
00658 }
00659 }
00660