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 "SundanceQuadratureIntegral.hpp"
00032 #include "SundanceGaussianQuadrature.hpp"
00033 #include "SundanceSpatialDerivSpecifier.hpp"
00034 #include "SundanceOut.hpp"
00035 #include "SundanceTabs.hpp"
00036 #include "Teuchos_TimeMonitor.hpp"
00037
00038 using namespace Sundance;
00039 using namespace Teuchos;
00040
00041 using std::endl;
00042 using std::setw;
00043 using std::setprecision;
00044
00045 extern "C"
00046 {
00047 int dgemm_(const char* transA, const char* transB,
00048 const int* M, const int *N, const int* K,
00049 const double* alpha,
00050 const double* A, const int* ldA,
00051 const double* B, const int* ldB,
00052 const double* beta,
00053 double* C, const int* ldC);
00054 }
00055
00056 static Time& quadratureTimer()
00057 {
00058 static RCP<Time> rtn
00059 = TimeMonitor::getNewTimer("quadrature");
00060 return *rtn;
00061 }
00062
00063
00064 QuadratureIntegral::QuadratureIntegral(int spatialDim,
00065 const CellType& maxCellType,
00066 int dim,
00067 const CellType& cellType,
00068 const QuadratureFamily& quad,
00069 bool isInternalBdry,
00070 const ParametrizedCurve& globalCurve,
00071 const Mesh& mesh,
00072 int verb)
00073 : QuadratureIntegralBase(spatialDim, maxCellType, dim, cellType, quad,
00074 isInternalBdry, globalCurve, mesh, verb),
00075 W_(),
00076 useSumFirstMethod_(true)
00077 {
00078 Tabs tab0(0);
00079
00080 SUNDANCE_MSG1(setupVerb(), tab0 << "QuadratureIntegral ctor for 0-form");
00081 if (setupVerb()) describe(Out::os());
00082
00083 SUNDANCE_MSG1(setupVerb(), tab0 << "quadrature family=" << quad);
00084
00085
00086 W_.resize(nFacetCases());
00087 quadPts_.resize(nFacetCases());
00088 quadWeights_.resize(nFacetCases());
00089 quad_ = quad;
00090
00091 for (int fc=0; fc<nFacetCases(); fc++)
00092 {
00093 Tabs tab2;
00094 SUNDANCE_MSG1(setupVerb(), tab2 << "facet case=" << fc
00095 << " of " << nFacetCases());
00096 Tabs tab3;
00097
00098 Array<double> quadWeights;
00099 Array<Point> quadPts;
00100 getQuad(quad, fc, quadPts, quadWeights);
00101 nQuad_ = quadPts.size();
00102 quadPts_[fc] = quadPts;
00103 quadWeights_[fc] = quadWeights;
00104
00105 W_[fc].resize(nQuad());
00106
00107 for (int q=0; q<nQuad(); q++)
00108 {
00109 W_[fc][q] = quadWeights[q];
00110 }
00111 }
00112 }
00113
00114
00115 QuadratureIntegral::QuadratureIntegral(int spatialDim,
00116 const CellType& maxCellType,
00117 int dim,
00118 const CellType& cellType,
00119 const BasisFamily& testBasis,
00120 int alpha,
00121 int testDerivOrder,
00122 const QuadratureFamily& quad,
00123 bool isInternalBdry,
00124 const ParametrizedCurve& globalCurve,
00125 const Mesh& mesh,
00126 int verb)
00127 : QuadratureIntegralBase(spatialDim, maxCellType, dim, cellType,
00128 testBasis, alpha, testDerivOrder, quad , isInternalBdry, globalCurve , mesh, verb),
00129 W_(),
00130 useSumFirstMethod_(true)
00131 {
00132 Tabs tab0;
00133
00134 SUNDANCE_MSG1(setupVerb(), tab0 << "QuadratureIntegral ctor for 1-form");
00135 if (setupVerb()) describe(Out::os());
00136 assertLinearForm();
00137
00138
00139 SUNDANCE_MSG1(setupVerb(), tab0 << "quadrature family=" << quad);
00140
00141 quad_ = quad;
00142 W_.resize(nFacetCases());
00143
00144 quadPts_.resize(nFacetCases());
00145 quadWeights_.resize(nFacetCases());
00146 W_ACI_F1_.resize(nFacetCases());
00147
00148 for (int fc=0; fc<nFacetCases(); fc++)
00149 {
00150 Tabs tab2;
00151
00152 SUNDANCE_MSG1(setupVerb(), tab2 << "facet case=" << fc
00153 << " of " << nFacetCases());
00154
00155
00156 Array<double> quadWeights;
00157 Array<Point> quadPts;
00158 getQuad(quad, fc, quadPts, quadWeights);
00159 nQuad_ = quadPts.size();
00160 quadPts_[fc] = quadPts;
00161 quadWeights_[fc] = quadWeights;
00162
00163 W_[fc].resize(nQuad() * nRefDerivTest() * nNodesTest());
00164
00165 SUNDANCE_MSG1(setupVerb(), tab2 << "num nodes for test function " << nNodesTest());
00166
00167 Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());
00168
00169 for (int r=0; r<nRefDerivTest(); r++)
00170 {
00171 Tabs tab3;
00172 SUNDANCE_MSG1(setupVerb(), tab3
00173 << "evaluating basis functions for ref deriv direction " << r);
00174 MultiIndex mi;
00175 testBasisVals[r].resize(testBasis.dim());
00176 if (testDerivOrder==1) mi[r] = 1;
00177 SpatialDerivSpecifier deriv(mi);
00178 testBasis.refEval(evalCellType(), quadPts, deriv,
00179 testBasisVals[r], setupVerb());
00180 }
00181
00182
00183
00184 int vecComp = 0;
00185 W_ACI_F1_[fc].resize(nQuad());
00186 for (int q=0; q<nQuad(); q++)
00187 {
00188 W_ACI_F1_[fc][q].resize(nRefDerivTest());
00189 for (int t=0; t<nRefDerivTest(); t++)
00190 {
00191 W_ACI_F1_[fc][q][t].resize(nNodesTest());
00192 for (int nt=0; nt<nNodesTest(); nt++)
00193 {
00194 wValue(fc, q, t, nt)
00195 = chop(quadWeights[q] * testBasisVals[t][vecComp][q][nt]) ;
00196 W_ACI_F1_[fc][q][t][nt] = chop(testBasisVals[t][vecComp][q][nt]);
00197 }
00198 }
00199 }
00200
00201 addFlops(2*nQuad()*nRefDerivTest()*nNodesTest());
00202 }
00203 }
00204
00205
00206
00207
00208 QuadratureIntegral::QuadratureIntegral(int spatialDim,
00209 const CellType& maxCellType,
00210 int dim,
00211 const CellType& cellType,
00212 const BasisFamily& testBasis,
00213 int alpha,
00214 int testDerivOrder,
00215 const BasisFamily& unkBasis,
00216 int beta,
00217 int unkDerivOrder,
00218 const QuadratureFamily& quad,
00219 bool isInternalBdry,
00220 const ParametrizedCurve& globalCurve,
00221 const Mesh& mesh,
00222 int verb)
00223 : QuadratureIntegralBase(spatialDim, maxCellType, dim, cellType,
00224 testBasis, alpha, testDerivOrder,
00225 unkBasis, beta, unkDerivOrder, quad , isInternalBdry, globalCurve , mesh , verb),
00226 W_(),
00227 useSumFirstMethod_(true)
00228 {
00229 Tabs tab0;
00230
00231 SUNDANCE_MSG1(setupVerb(), tab0 << "QuadratureIntegral ctor for 2-form");
00232 if (setupVerb()) describe(Out::os());
00233 assertBilinearForm();
00234 quad_ = quad;
00235
00236 W_.resize(nFacetCases());
00237 W_ACI_F2_.resize(nFacetCases());
00238
00239 quadPts_.resize(nFacetCases());
00240 quadWeights_.resize(nFacetCases());
00241
00242 for (int fc=0; fc<nFacetCases(); fc++)
00243 {
00244
00245 Array<double> quadWeights;
00246 Array<Point> quadPts;
00247 getQuad(quad, fc, quadPts, quadWeights);
00248
00249 quadPts_[fc] = quadPts;
00250 quadWeights_[fc] = quadWeights;
00251 nQuad_ = quadPts.size();
00252
00253 W_[fc].resize(nQuad() * nRefDerivTest() * nNodesTest()
00254 * nRefDerivUnk() * nNodesUnk());
00255
00256
00257
00258 Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());
00259 Array<Array<Array<Array<double> > > > unkBasisVals(nRefDerivUnk());
00260
00261
00262 for (int r=0; r<nRefDerivTest(); r++)
00263 {
00264 testBasisVals[r].resize(testBasis.dim());
00265 MultiIndex mi;
00266 if (testDerivOrder==1) mi[r] = 1;
00267 SpatialDerivSpecifier deriv(mi);
00268 testBasis.refEval(evalCellType(), quadPts, deriv,
00269 testBasisVals[r], setupVerb());
00270 }
00271
00272 for (int r=0; r<nRefDerivUnk(); r++)
00273 {
00274 unkBasisVals[r].resize(unkBasis.dim());
00275 MultiIndex mi;
00276 if (unkDerivOrder==1) mi[r] = 1;
00277 SpatialDerivSpecifier deriv(mi);
00278 unkBasis.refEval(evalCellType(),
00279 quadPts, deriv, unkBasisVals[r], setupVerb());
00280 }
00281
00282
00283 int vecComp = 0;
00284
00285 W_ACI_F2_[fc].resize(nQuad());
00286 for (int q=0; q<nQuad(); q++)
00287 {
00288 W_ACI_F2_[fc][q].resize(nRefDerivTest());
00289 for (int t=0; t<nRefDerivTest(); t++)
00290 {
00291 W_ACI_F2_[fc][q][t].resize(nNodesTest());
00292 for (int nt=0; nt<nNodesTest(); nt++)
00293 {
00294 W_ACI_F2_[fc][q][t][nt].resize(nRefDerivUnk());
00295 for (int u=0; u<nRefDerivUnk(); u++)
00296 {
00297 W_ACI_F2_[fc][q][t][nt][u].resize(nNodesUnk());
00298 for (int nu=0; nu<nNodesUnk(); nu++)
00299 {
00300 wValue(fc, q, t, nt, u, nu)
00301 = chop(quadWeights[q] * testBasisVals[t][vecComp][q][nt]
00302 * unkBasisVals[u][vecComp][q][nu]);
00303 W_ACI_F2_[fc][q][t][nt][u][nu] =
00304 chop(testBasisVals[t][vecComp][q][nt] * unkBasisVals[u][vecComp][q][nu]);
00305 }
00306 }
00307 }
00308 }
00309 }
00310
00311 addFlops(3*nQuad()*nRefDerivTest()*nNodesTest()*nRefDerivUnk()*nNodesUnk()
00312 + W_[fc].size());
00313 for (int i=0; i<W_[fc].size(); i++) W_[fc][i] = chop(W_[fc][i]);
00314
00315 }
00316
00317 }
00318
00319
00320 void QuadratureIntegral::transformZeroForm(const CellJacobianBatch& JTrans,
00321 const CellJacobianBatch& JVol,
00322 const Array<int>& isLocalFlag,
00323 const Array<int>& facetIndex,
00324 const RCP<Array<int> >& cellLIDs,
00325 const double* const coeff,
00326 RCP<Array<double> >& A) const
00327 {
00328 TimeMonitor timer(quadratureTimer());
00329 Tabs tabs;
00330 TEST_FOR_EXCEPTION(order() != 0, InternalError,
00331 "QuadratureIntegral::transformZeroForm() called "
00332 "for form of order " << order());
00333
00334 TEST_FOR_EXCEPTION( (int) isLocalFlag.size() != 0
00335 && (int) isLocalFlag.size() != JVol.numCells(),
00336 RuntimeError,
00337 "mismatch between isLocalFlag.size()=" << isLocalFlag.size()
00338 << " and JVol.numCells()=" << JVol.numCells());
00339
00340 bool checkLocalFlag = (int) isLocalFlag.size() != 0;
00341
00342 const Array<int>* cellLID = cellLIDs.get();
00343
00344 SUNDANCE_MSG1(integrationVerb(), tabs << "doing zero form by quadrature");
00345 double& a = (*A)[0];
00346 SUNDANCE_MSG5(integrationVerb(), tabs << "input A=");
00347 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00348 double* coeffPtr = (double*) coeff;
00349
00350 int flops = 0;
00351 if (nFacetCases()==1)
00352 {
00353 if (globalCurve().isCurveValid())
00354 {
00355
00356 Array<double> quadWeightsTmp = quadWeights_[0];
00357 Array<Point> quadPointsTmp = quadPts_[0];
00358 bool isCutByCurve;
00359
00360 const Array<double>& w = W_[0];
00361 for (int c=0; c<JVol.numCells(); c++)
00362 {
00363 if (checkLocalFlag && !isLocalFlag[c])
00364 {
00365 coeffPtr += nQuad();
00366 continue;
00367 }
00368 double detJ = fabs(JVol.detJ()[c]);
00369 int fc = 0;
00370
00371 quadWeightsTmp = quadWeights_[fc];
00372 quadPointsTmp = quadPts_[fc];
00373 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c], fc ,mesh(),
00374 globalCurve(), quadPointsTmp, quadWeightsTmp, isCutByCurve);
00375
00376 if (isCutByCurve)
00377 {
00378 for (int q=0; q<nQuad(); q++, coeffPtr++)
00379 {
00380 a += quadWeightsTmp[q]*(*coeffPtr)*detJ;
00381 }
00382 flops += 3*nQuad();
00383 }
00384 else
00385 {
00386 for (int q=0; q<nQuad(); q++, coeffPtr++)
00387 {
00388 a += w[q]*(*coeffPtr)*detJ;
00389 }
00390 flops += 3*nQuad();
00391 }
00392
00393 }
00394 }
00395 else
00396 {
00397 const Array<double>& w = W_[0];
00398 for (int c=0; c<JVol.numCells(); c++)
00399 {
00400 if (checkLocalFlag && !isLocalFlag[c])
00401 {
00402 coeffPtr += nQuad();
00403 continue;
00404 }
00405 double detJ = fabs(JVol.detJ()[c]);
00406
00407 for (int q=0; q<nQuad(); q++, coeffPtr++)
00408 {
00409 a += w[q]*(*coeffPtr)*detJ;
00410 }
00411 flops += 3*nQuad();
00412 }
00413 }
00414 }
00415 else
00416 {
00417 if (globalCurve().isCurveValid())
00418 {
00419 Array<double> quadWeightsTmp = quadWeights_[0];
00420 Array<Point> quadPointsTmp = quadPts_[0];
00421 bool isCutByCurve;
00422
00423 for (int c=0; c<JVol.numCells(); c++)
00424 {
00425 if (checkLocalFlag && !isLocalFlag[c])
00426 {
00427 coeffPtr += nQuad();
00428 continue;
00429 }
00430 double detJ = fabs(JVol.detJ()[c]);
00431 int fc = facetIndex[c];
00432 const Array<double>& w = W_[facetIndex[c]];
00433
00434
00435 quadWeightsTmp = quadWeights_[fc];
00436 quadPointsTmp = quadPts_[fc];
00437 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c], fc ,mesh(),
00438 globalCurve(), quadPointsTmp, quadWeightsTmp, isCutByCurve);
00439
00440 if (isCutByCurve){
00441 for (int q=0; q<nQuad(); q++, coeffPtr++)
00442 a += quadWeightsTmp[q]*(*coeffPtr)*detJ;
00443 flops += 3*nQuad();
00444 }
00445 else{
00446 for (int q=0; q<nQuad(); q++, coeffPtr++)
00447 a += w[q]*(*coeffPtr)*detJ;
00448 flops += 3*nQuad();
00449 }
00450 }
00451 }
00452 else
00453 {
00454 for (int c=0; c<JVol.numCells(); c++)
00455 {
00456 if (checkLocalFlag && !isLocalFlag[c])
00457 {
00458 coeffPtr += nQuad();
00459 continue;
00460 }
00461 double detJ = fabs(JVol.detJ()[c]);
00462 const Array<double>& w = W_[facetIndex[c]];
00463
00464 for (int q=0; q<nQuad(); q++, coeffPtr++)
00465 {
00466 a += w[q]*(*coeffPtr)*detJ;
00467 }
00468 flops += 3*nQuad();
00469 }
00470 }
00471 }
00472 SUNDANCE_MSG5(integrationVerb(), tabs << "output A = ");
00473 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00474
00475 addFlops(flops);
00476 }
00477
00478
00479 void QuadratureIntegral::transformOneForm(const CellJacobianBatch& JTrans,
00480 const CellJacobianBatch& JVol,
00481 const Array<int>& facetIndex,
00482 const RCP<Array<int> >& cellLIDs,
00483 const double* const coeff,
00484 RCP<Array<double> >& A) const
00485 {
00486 TimeMonitor timer(quadratureTimer());
00487 Tabs tabs;
00488 TEST_FOR_EXCEPTION(order() != 1, InternalError,
00489 "QuadratureIntegral::transformOneForm() called for form "
00490 "of order " << order());
00491 SUNDANCE_MSG2(integrationVerb(), tabs << "doing one form by quadrature");
00492 int flops = 0;
00493
00494 const Array<int>* cellLID = cellLIDs.get();
00495
00496
00497
00498
00499 if (testDerivOrder() == 0)
00500 {
00501 double* aPtr = &((*A)[0]);
00502 SUNDANCE_MSG5(integrationVerb(), tabs << "input A = ");
00503 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00504
00505 double* coeffPtr = (double*) coeff;
00506 int offset = 0 ;
00507
00508 if (globalCurve().isCurveValid())
00509 {
00510
00511 Array<double> quadWeightsTmp = quadWeights_[0];
00512 Array<Point> quadPointsTmp = quadPts_[0];
00513 bool isCutByCurve;
00514
00515 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00516 {
00517 Tabs tab2;
00518 double detJ = fabs(JVol.detJ()[c]);
00519 int fc = 0;
00520 if (nFacetCases() != 1) fc = facetIndex[c];
00521 const Array<double>& w = W_[fc];
00522 SUNDANCE_MSG4(integrationVerb(), tab2 << "c=" << c << " detJ=" << detJ);
00523
00524 quadWeightsTmp = quadWeights_[fc];
00525 quadPointsTmp = quadPts_[fc];
00526
00527 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00528 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00529 if (isCutByCurve){
00530 Array<double> wi;
00531 wi.resize(nQuad() * nNodes());
00532 for (int ii = 0 ; ii < wi.size() ; ii++ ) { wi[ii] = 0.0; }
00533 for (int n = 0 ; n < nNodes() ; n++)
00534 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00535
00536 wi[n + nNodes()*q] +=
00537 chop(quadWeightsTmp[q] * W_ACI_F1_[fc][q][0][n]);
00538
00539 for (int q=0; q<nQuad(); q++, coeffPtr++){
00540 double f = (*coeffPtr)*detJ;
00541 for (int n=0; n<nNodes(); n++){
00542 aPtr[offset+n] += f*wi[n + nNodes()*q];
00543 }
00544 }
00545 }
00546 else
00547 {
00548 for (int q=0; q<nQuad(); q++, coeffPtr++){
00549 double f = (*coeffPtr)*detJ;
00550 for (int n=0; n<nNodes(); n++){
00551 aPtr[offset+n] += f*w[n + nNodes()*q];
00552 }
00553 }
00554 }
00555 }
00556 }
00557 else
00558 {
00559 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00560 {
00561 Tabs tab2;
00562 double detJ = fabs(JVol.detJ()[c]);
00563 int fc = 0;
00564 if (nFacetCases() != 1) fc = facetIndex[c];
00565 const Array<double>& w = W_[fc];
00566 SUNDANCE_MSG4(integrationVerb(), tab2 << "c=" << c << " detJ=" << detJ);
00567
00568 for (int q=0; q<nQuad(); q++, coeffPtr++)
00569 {
00570 Tabs tab3;
00571 double f = (*coeffPtr)*detJ;
00572 SUNDANCE_MSG4(integrationVerb(), tab3 << "q=" << q << " coeff=" <<
00573 *coeffPtr << " coeff*detJ=" << f);
00574 for (int n=0; n<nNodes(); n++)
00575 {
00576 Tabs tab4;
00577 SUNDANCE_MSG4(integrationVerb(), tab4 << "n=" << n << " w=" <<
00578 w[n + nNodes()*q]);
00579 aPtr[offset+n] += f*w[n + nNodes()*q];
00580 }
00581 }
00582 }
00583 }
00584 if (integrationVerb() >= 4)
00585 {
00586 Tabs tab2;
00587 Out::os() << tab2 << "integration results on cell:" << std::endl;
00588 Out::os() << tab2 << setw(10) << "n" << setw(12) << "I_n" << std::endl;
00589 for (int n=0; n<nNodes(); n++)
00590 {
00591 Out::os() << tab2 << setw(10) << n
00592 << setw(12) << setprecision(6) << aPtr[offset+n] << std::endl;
00593 }
00594 }
00595
00596 SUNDANCE_MSG5(integrationVerb(), tabs << "output A = ");
00597 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00598 addFlops( JVol.numCells() * (1 + nQuad() * (1 + 2*nNodes())) );
00599 }
00600 else
00601 {
00602
00603
00604
00605
00606 createOneFormTransformationMatrix(JTrans, JVol);
00607
00608 SUNDANCE_MSG4(transformVerb(),
00609 Tabs() << "transformation matrix=" << G(alpha()));
00610
00611 double* GPtr = &(G(alpha())[0]);
00612
00613 if (useSumFirstMethod())
00614 {
00615 transformSummingFirst(JVol.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00616 }
00617 else
00618 {
00619 transformSummingLast(JVol.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00620 }
00621 }
00622 addFlops(flops);
00623 }
00624
00625
00626 void QuadratureIntegral::transformTwoForm(const CellJacobianBatch& JTrans,
00627 const CellJacobianBatch& JVol,
00628 const Array<int>& facetIndex,
00629 const RCP<Array<int> >& cellLIDs,
00630 const double* const coeff,
00631 RCP<Array<double> >& A) const
00632 {
00633 TimeMonitor timer(quadratureTimer());
00634 Tabs tabs;
00635 TEST_FOR_EXCEPTION(order() != 2, InternalError,
00636 "QuadratureIntegral::transformTwoForm() called for form "
00637 "of order " << order());
00638 SUNDANCE_MSG2(integrationVerb(), tabs << "doing two form by quadrature");
00639
00640 const Array<int>* cellLID = cellLIDs.get();
00641
00642
00643
00644
00645 if (testDerivOrder() == 0 && unkDerivOrder() == 0)
00646 {
00647 double* aPtr = &((*A)[0]);
00648 double* coeffPtr = (double*) coeff;
00649 int offset = 0 ;
00650
00651 if (globalCurve().isCurveValid())
00652 {
00653
00654 Array<double> quadWeightsTmp = quadWeights_[0];
00655 Array<Point> quadPointsTmp = quadPts_[0];
00656 bool isCutByCurve;
00657
00658 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00659 {
00660 double detJ = fabs(JVol.detJ()[c]);
00661 int fc = 0;
00662 if (nFacetCases() != 1) fc = facetIndex[c];
00663 const Array<double>& w = W_[fc];
00664
00665 quadWeightsTmp = quadWeights_[fc];
00666 quadPointsTmp = quadPts_[fc];
00667
00668 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00669 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00670 if (isCutByCurve){
00671 Array<double> wi;
00672 wi.resize(nQuad() * nNodesTest() *nNodesUnk() );
00673 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00674 for (int nt = 0 ; nt < nNodesTest() ; nt++){
00675 for (int nu=0; nu<nNodesUnk(); nu++)
00676 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00677
00678 wi[nu + nNodesUnk()*(nt + nNodesTest()*q)] +=
00679 chop(quadWeightsTmp[q] * W_ACI_F2_[fc][q][0][nt][0][nu]);
00680 }
00681 for (int q=0; q<nQuad(); q++, coeffPtr++){
00682 double f = (*coeffPtr)*detJ;
00683 for (int n=0; n<nNodes(); n++){
00684 aPtr[offset+n] += f*wi[n + nNodes()*q];
00685 }
00686 }
00687 }
00688 else
00689 {
00690 for (int q=0; q<nQuad(); q++, coeffPtr++){
00691 double f = (*coeffPtr)*detJ;
00692 for (int n=0; n<nNodes(); n++){
00693 aPtr[offset+n] += f*w[n + nNodes()*q];
00694 }
00695 }
00696 }
00697 }
00698 }
00699 else
00700 {
00701
00702 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00703 {
00704 double detJ = fabs(JVol.detJ()[c]);
00705 int fc = 0;
00706 if (nFacetCases() != 1) fc = facetIndex[c];
00707 const Array<double>& w = W_[fc];
00708 for (int q=0; q<nQuad(); q++, coeffPtr++)
00709 {
00710 double f = (*coeffPtr)*detJ;
00711 for (int n=0; n<nNodes(); n++)
00712 {
00713 aPtr[offset+n] += f*w[n + nNodes()*q];
00714 }
00715 }
00716 }
00717 }
00718 addFlops( JVol.numCells() * (1 + nQuad() * (1 + 2*nNodes())) );
00719 }
00720 else
00721 {
00722 createTwoFormTransformationMatrix(JTrans, JVol);
00723 double* GPtr;
00724
00725 if (testDerivOrder() == 0)
00726 {
00727 GPtr = &(G(beta())[0]);
00728 SUNDANCE_MSG2(transformVerb(),
00729 Tabs() << "transformation matrix=" << G(beta()));
00730 }
00731 else if (unkDerivOrder() == 0)
00732 {
00733 GPtr = &(G(alpha())[0]);
00734 SUNDANCE_MSG2(transformVerb(),
00735 Tabs() << "transformation matrix=" << G(alpha()));
00736 }
00737 else
00738 {
00739 GPtr = &(G(alpha(), beta())[0]);
00740 SUNDANCE_MSG2(transformVerb(),
00741 Tabs() << "transformation matrix="
00742 << G(alpha(),beta()));
00743 }
00744
00745
00746 if (useSumFirstMethod())
00747 {
00748 transformSummingFirst(JTrans.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00749 }
00750 else
00751 {
00752 transformSummingLast(JTrans.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00753 }
00754 }
00755 }
00756
00757 void QuadratureIntegral
00758 ::transformSummingFirst(int nCells,
00759 const Array<int>& facetIndex,
00760 const RCP<Array<int> >& cellLIDs,
00761 const double* const GPtr,
00762 const double* const coeff,
00763 RCP<Array<double> >& A) const
00764 {
00765 double* aPtr = &((*A)[0]);
00766 double* coeffPtr = (double*) coeff;
00767 const Array<int>* cellLID = cellLIDs.get();
00768
00769 int transSize = 0;
00770 if (order()==2)
00771 {
00772 transSize = nRefDerivTest() * nRefDerivUnk();
00773 }
00774 else
00775 {
00776 transSize = nRefDerivTest();
00777 }
00778
00779
00780 static Array<double> sumWorkspace;
00781
00782 int swSize = transSize * nNodes();
00783 sumWorkspace.resize(swSize);
00784
00785
00786
00787
00788
00789
00790
00791
00792
00793
00794
00795 if (globalCurve().isCurveValid())
00796 {
00797
00798 Array<double> quadWeightsTmp = quadWeights_[0];
00799 Array<Point> quadPointsTmp = quadPts_[0];
00800 bool isCutByCurve;
00801
00802 for (int c=0; c<nCells; c++)
00803 {
00804
00805 for (int i=0; i<swSize; i++) sumWorkspace[i]=0.0;
00806 int fc = 0;
00807 if (nFacetCases() > 1) fc = facetIndex[c];
00808
00809 const Array<double>& w = W_[fc];
00810
00811 quadWeightsTmp = quadWeights_[fc];
00812 quadPointsTmp = quadPts_[fc];
00813
00814 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00815 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00816 if (isCutByCurve){
00817 Array<double> wi;
00818 if (order()==1){
00819 wi.resize(nQuad() * nRefDerivTest() * nNodesTest());
00820 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00821 for (int n = 0 ; n < nNodes() ; n++)
00822 for (int t=0; t<nRefDerivTest(); t++)
00823 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00824
00825 wi[n + nNodes()*(t + nRefDerivTest()*q)] +=
00826 chop(quadWeightsTmp[q] * W_ACI_F1_[fc][q][t][n]);
00827 }else{
00828 wi.resize(nQuad() * nRefDerivTest() * nNodesTest()
00829 * nRefDerivUnk() * nNodesUnk() );
00830 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00831
00832 for (int t=0; t<nRefDerivTest(); t++)
00833 for (int nt = 0 ; nt < nNodesTest() ; nt++){
00834 for (int u=0; u<nRefDerivUnk(); u++)
00835 for (int nu=0; nu<nNodesUnk(); nu++)
00836 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00837
00838
00839 wi[nu + nNodesUnk()*(nt + nNodesTest()*(u + nRefDerivUnk()*(t + nRefDerivTest()*q)))] +=
00840 chop(quadWeightsTmp[q] * W_ACI_F2_[fc][q][t][nt][u][nu]);
00841 }
00842 }
00843 for (int q=0; q<nQuad(); q++, coeffPtr++){
00844 double f = (*coeffPtr);
00845 for (int n=0; n<swSize; n++){
00846 sumWorkspace[n] += f*wi[n + q*swSize];
00847 }
00848 }
00849 }
00850 else{
00851 for (int q=0; q<nQuad(); q++, coeffPtr++)
00852 {
00853 double f = (*coeffPtr);
00854 for (int n=0; n<swSize; n++)
00855 {
00856 sumWorkspace[n] += f*w[n + q*swSize];
00857 }
00858 }
00859 }
00860 }
00861 }
00862 else
00863 {
00864 for (int c=0; c<nCells; c++)
00865 {
00866
00867 for (int i=0; i<swSize; i++) sumWorkspace[i]=0.0;
00868 int fc = 0;
00869 if (nFacetCases() > 1) fc = facetIndex[c];
00870
00871 const Array<double>& w = W_[fc];
00872
00873 for (int q=0; q<nQuad(); q++, coeffPtr++)
00874 {
00875 double f = (*coeffPtr);
00876 for (int n=0; n<swSize; n++)
00877 {
00878 sumWorkspace[n] += f*w[n + q*swSize];
00879 }
00880 }
00881
00882
00883 const double * const gCell = &(GPtr[transSize*c]);
00884 double* aCell = aPtr + nNodes()*c;
00885 for (int i=0; i<nNodes(); i++)
00886 {
00887 for (int j=0; j<transSize; j++)
00888 {
00889 aCell[i] += sumWorkspace[nNodes()*j + i] * gCell[j];
00890 }
00891 }
00892 }
00893 }
00894 int flops = 2*(nCells * nNodes() * transSize) * (nQuad() + 1) ;
00895 addFlops(flops);
00896 }
00897
00898 void QuadratureIntegral
00899 ::transformSummingLast(int nCells,
00900 const Array<int>& facetIndex,
00901 const RCP<Array<int> >& cellLIDs,
00902 const double* const GPtr,
00903 const double* const coeff,
00904 RCP<Array<double> >& A) const
00905 {
00906 double* aPtr = &((*A)[0]);
00907 int transSize = 0;
00908 const Array<int>* cellLID = cellLIDs.get();
00909
00910 if (order()==2)
00911 {
00912 transSize = nRefDerivTest() * nRefDerivUnk();
00913 }
00914 else
00915 {
00916 transSize = nRefDerivTest();
00917 }
00918
00919
00920
00921 static Array<double> jWorkspace;
00922 jWorkspace.resize(transSize);
00923
00924
00925
00926
00927
00928
00929
00930
00931
00932
00933 if (globalCurve().isCurveValid()){
00934
00935 Array<double> quadWeightsTmp = quadWeights_[0];
00936 Array<Point> quadPointsTmp = quadPts_[0];
00937 bool isCutByCurve;
00938
00939 for (int c=0; c<nCells; c++)
00940 {
00941 const double* const gCell = &(GPtr[transSize*c]);
00942 double* aCell = aPtr + nNodes()*c;
00943 int fc = 0;
00944 if (nFacetCases() > 1) fc = facetIndex[c];
00945 const Array<double>& w = W_[fc];
00946
00947 quadWeightsTmp = quadWeights_[fc];
00948 quadPointsTmp = quadPts_[fc];
00949
00950 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00951 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00952 if (isCutByCurve){
00953 Array<double> wi;
00954 if (order()==1){
00955 wi.resize(nQuad() * nRefDerivTest() * nNodesTest());
00956 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00957 for (int n = 0 ; n < nNodes() ; n++)
00958 for (int t=0; t<nRefDerivTest(); t++)
00959 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00960
00961 wi[n + nNodes()*(t + nRefDerivTest()*q)] +=
00962 chop(quadWeightsTmp[q] * W_ACI_F1_[fc][q][t][n]);
00963 }else{
00964 wi.resize(nQuad() * nRefDerivTest() * nNodesTest()
00965 * nRefDerivUnk() * nNodesUnk() );
00966 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00967
00968 for (int t=0; t<nRefDerivTest(); t++)
00969 for (int nt = 0 ; nt < nNodesTest() ; nt++){
00970 for (int u=0; u<nRefDerivUnk(); u++)
00971 for (int nu=0; nu<nNodesUnk(); nu++)
00972 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00973
00974
00975 wi[nu + nNodesUnk()*(nt + nNodesTest()*(u + nRefDerivUnk()*(t + nRefDerivTest()*q)))] +=
00976 chop(quadWeightsTmp[q] * W_ACI_F2_[fc][q][t][nt][u][nu]);
00977 }
00978 }
00979 for (int q=0; q<nQuad(); q++){
00980 double f = coeff[c*nQuad() + q];
00981 for (int t=0; t<transSize; t++) jWorkspace[t]=f*gCell[t];
00982
00983 for (int n=0; n<nNodes(); n++){
00984 for (int t=0; t<transSize; t++){
00985 aCell[n] += jWorkspace[t]*wi[n + nNodes()*(t + transSize*q)];
00986 }
00987 }
00988 }
00989 }
00990 else{
00991 for (int q=0; q<nQuad(); q++){
00992 double f = coeff[c*nQuad() + q];
00993 for (int t=0; t<transSize; t++) jWorkspace[t]=f*gCell[t];
00994
00995 for (int n=0; n<nNodes(); n++){
00996 for (int t=0; t<transSize; t++){
00997 aCell[n] += jWorkspace[t]*w[n + nNodes()*(t + transSize*q)];
00998 }
00999 }
01000 }
01001 }
01002 }
01003 }
01004 else
01005 {
01006 for (int c=0; c<nCells; c++)
01007 {
01008 const double* const gCell = &(GPtr[transSize*c]);
01009 double* aCell = aPtr + nNodes()*c;
01010 int fc = 0;
01011 if (nFacetCases() > 1) fc = facetIndex[c];
01012 const Array<double>& w = W_[fc];
01013 for (int q=0; q<nQuad(); q++)
01014 {
01015 double f = coeff[c*nQuad() + q];
01016 for (int t=0; t<transSize; t++) jWorkspace[t]=f*gCell[t];
01017
01018 for (int n=0; n<nNodes(); n++)
01019 {
01020 for (int t=0; t<transSize; t++)
01021 {
01022 aCell[n] += jWorkspace[t]*w[n + nNodes()*(t + transSize*q)];
01023 }
01024 }
01025 }
01026 }
01027 }
01028 int flops = nCells * nQuad() * transSize * (1 + 2*nNodes()) ;
01029 addFlops(flops);
01030 }