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