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