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 "PlayaTabs.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 TEUCHOS_TEST_FOR_EXCEPTION(order() != 0, std::logic_error,
00331 "QuadratureIntegral::transformZeroForm() called "
00332 "for form of order " << order());
00333
00334 TEUCHOS_TEST_FOR_EXCEPTION( (int) isLocalFlag.size() != 0
00335 && (int) isLocalFlag.size() != JVol.numCells(),
00336 std::runtime_error,
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 Tabs tab1;
00457 SUNDANCE_MSG4(integrationVerb(), tab1 << "cell #" << c
00458 << " of " << JVol.numCells());
00459 if (checkLocalFlag && !isLocalFlag[c])
00460 {
00461 Tabs tab2;
00462 SUNDANCE_MSG4(integrationVerb(), tab2
00463 << "skipped -- is not local cell");
00464 coeffPtr += nQuad();
00465 continue;
00466 }
00467 double detJ = fabs(JVol.detJ()[c]);
00468 const Array<double>& w = W_[facetIndex[c]];
00469
00470 for (int q=0; q<nQuad(); q++, coeffPtr++)
00471 {
00472 Tabs tab3;
00473 SUNDANCE_MSG4(integrationVerb(),
00474 tab3 << "q=" << q << "\t w=" << w[q]
00475 << "\t\t coeff=" << *coeffPtr << "\t\t |J|=" << detJ);
00476 a += w[q]*(*coeffPtr)*detJ;
00477 }
00478 flops += 3*nQuad();
00479 }
00480 }
00481 }
00482 SUNDANCE_MSG5(integrationVerb(), tabs << "output A = ");
00483 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00484
00485 addFlops(flops);
00486 }
00487
00488
00489 void QuadratureIntegral::transformOneForm(const CellJacobianBatch& JTrans,
00490 const CellJacobianBatch& JVol,
00491 const Array<int>& facetIndex,
00492 const RCP<Array<int> >& cellLIDs,
00493 const double* const coeff,
00494 RCP<Array<double> >& A) const
00495 {
00496 TimeMonitor timer(quadratureTimer());
00497 Tabs tabs;
00498 TEUCHOS_TEST_FOR_EXCEPTION(order() != 1, std::logic_error,
00499 "QuadratureIntegral::transformOneForm() called for form "
00500 "of order " << order());
00501 SUNDANCE_MSG2(integrationVerb(), tabs << "doing one form by quadrature");
00502 int flops = 0;
00503
00504 const Array<int>* cellLID = cellLIDs.get();
00505
00506
00507
00508
00509 if (testDerivOrder() == 0)
00510 {
00511 double* aPtr = &((*A)[0]);
00512 SUNDANCE_MSG5(integrationVerb(), tabs << "input A = ");
00513 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00514
00515 double* coeffPtr = (double*) coeff;
00516 int offset = 0 ;
00517
00518 if (globalCurve().isCurveValid())
00519 {
00520
00521 Array<double> quadWeightsTmp = quadWeights_[0];
00522 Array<Point> quadPointsTmp = quadPts_[0];
00523 bool isCutByCurve;
00524
00525 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00526 {
00527 Tabs tab2;
00528 double detJ = fabs(JVol.detJ()[c]);
00529 int fc = 0;
00530 if (nFacetCases() != 1) fc = facetIndex[c];
00531 const Array<double>& w = W_[fc];
00532 SUNDANCE_MSG4(integrationVerb(), tab2 << "c=" << c << " detJ=" << detJ);
00533
00534 quadWeightsTmp = quadWeights_[fc];
00535 quadPointsTmp = quadPts_[fc];
00536
00537 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00538 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00539 if (isCutByCurve){
00540 Array<double> wi;
00541 wi.resize(nQuad() * nNodes());
00542 for (int ii = 0 ; ii < wi.size() ; ii++ ) { wi[ii] = 0.0; }
00543 for (int n = 0 ; n < nNodes() ; n++)
00544 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00545
00546 wi[n + nNodes()*q] +=
00547 chop(quadWeightsTmp[q] * W_ACI_F1_[fc][q][0][n]);
00548
00549 for (int q=0; q<nQuad(); q++, coeffPtr++){
00550 double f = (*coeffPtr)*detJ;
00551 for (int n=0; n<nNodes(); n++){
00552 aPtr[offset+n] += f*wi[n + nNodes()*q];
00553 }
00554 }
00555 }
00556 else
00557 {
00558 for (int q=0; q<nQuad(); q++, coeffPtr++){
00559 double f = (*coeffPtr)*detJ;
00560 for (int n=0; n<nNodes(); n++){
00561 aPtr[offset+n] += f*w[n + nNodes()*q];
00562 }
00563 }
00564 }
00565 }
00566 }
00567 else
00568 {
00569 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00570 {
00571 Tabs tab2;
00572 double detJ = fabs(JVol.detJ()[c]);
00573 int fc = 0;
00574 if (nFacetCases() != 1) fc = facetIndex[c];
00575 const Array<double>& w = W_[fc];
00576 SUNDANCE_MSG4(integrationVerb(), tab2 << "c=" << c << " detJ=" << detJ);
00577
00578 for (int q=0; q<nQuad(); q++, coeffPtr++)
00579 {
00580 Tabs tab3;
00581 double f = (*coeffPtr)*detJ;
00582 SUNDANCE_MSG4(integrationVerb(), tab3 << "q=" << q << " coeff=" <<
00583 *coeffPtr << " coeff*detJ=" << f);
00584 for (int n=0; n<nNodes(); n++)
00585 {
00586 Tabs tab4;
00587 SUNDANCE_MSG4(integrationVerb(), tab4 << "n=" << n << " w=" <<
00588 w[n + nNodes()*q]);
00589 aPtr[offset+n] += f*w[n + nNodes()*q];
00590 }
00591 }
00592 }
00593 }
00594 if (integrationVerb() >= 4)
00595 {
00596 Tabs tab2;
00597 Out::os() << tab2 << "integration results on cell:" << std::endl;
00598 Out::os() << tab2 << setw(10) << "n" << setw(12) << "I_n" << std::endl;
00599 for (int n=0; n<nNodes(); n++)
00600 {
00601 Out::os() << tab2 << setw(10) << n
00602 << setw(12) << setprecision(6) << aPtr[offset+n] << std::endl;
00603 }
00604 }
00605
00606 SUNDANCE_MSG5(integrationVerb(), tabs << "output A = ");
00607 if (integrationVerb() >= 5) writeTable(Out::os(), tabs, *A, 6);
00608 addFlops( JVol.numCells() * (1 + nQuad() * (1 + 2*nNodes())) );
00609 }
00610 else
00611 {
00612
00613
00614
00615
00616 createOneFormTransformationMatrix(JTrans, JVol);
00617
00618 SUNDANCE_MSG4(transformVerb(),
00619 Tabs() << "transformation matrix=" << G(alpha()));
00620
00621 double* GPtr = &(G(alpha())[0]);
00622
00623 if (useSumFirstMethod())
00624 {
00625 transformSummingFirst(JVol.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00626 }
00627 else
00628 {
00629 transformSummingLast(JVol.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00630 }
00631 }
00632 addFlops(flops);
00633 }
00634
00635
00636 void QuadratureIntegral::transformTwoForm(const CellJacobianBatch& JTrans,
00637 const CellJacobianBatch& JVol,
00638 const Array<int>& facetIndex,
00639 const RCP<Array<int> >& cellLIDs,
00640 const double* const coeff,
00641 RCP<Array<double> >& A) const
00642 {
00643 TimeMonitor timer(quadratureTimer());
00644 Tabs tabs;
00645 TEUCHOS_TEST_FOR_EXCEPTION(order() != 2, std::logic_error,
00646 "QuadratureIntegral::transformTwoForm() called for form "
00647 "of order " << order());
00648 SUNDANCE_MSG2(integrationVerb(), tabs << "doing two form by quadrature");
00649
00650 const Array<int>* cellLID = cellLIDs.get();
00651
00652
00653
00654
00655 if (testDerivOrder() == 0 && unkDerivOrder() == 0)
00656 {
00657 double* aPtr = &((*A)[0]);
00658 double* coeffPtr = (double*) coeff;
00659 int offset = 0 ;
00660
00661 if (globalCurve().isCurveValid())
00662 {
00663
00664 Array<double> quadWeightsTmp = quadWeights_[0];
00665 Array<Point> quadPointsTmp = quadPts_[0];
00666 bool isCutByCurve;
00667
00668 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00669 {
00670 double detJ = fabs(JVol.detJ()[c]);
00671 int fc = 0;
00672 if (nFacetCases() != 1) fc = facetIndex[c];
00673 const Array<double>& w = W_[fc];
00674
00675 quadWeightsTmp = quadWeights_[fc];
00676 quadPointsTmp = quadPts_[fc];
00677
00678 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00679 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00680 if (isCutByCurve){
00681 Array<double> wi;
00682 wi.resize(nQuad() * nNodesTest() *nNodesUnk() );
00683 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00684 for (int nt = 0 ; nt < nNodesTest() ; nt++){
00685 for (int nu=0; nu<nNodesUnk(); nu++)
00686 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00687
00688 wi[nu + nNodesUnk()*(nt + nNodesTest()*q)] +=
00689 chop(quadWeightsTmp[q] * W_ACI_F2_[fc][q][0][nt][0][nu]);
00690 }
00691 for (int q=0; q<nQuad(); q++, coeffPtr++){
00692 double f = (*coeffPtr)*detJ;
00693 for (int n=0; n<nNodes(); n++){
00694 aPtr[offset+n] += f*wi[n + nNodes()*q];
00695 }
00696 }
00697 }
00698 else
00699 {
00700 for (int q=0; q<nQuad(); q++, coeffPtr++){
00701 double f = (*coeffPtr)*detJ;
00702 for (int n=0; n<nNodes(); n++){
00703 aPtr[offset+n] += f*w[n + nNodes()*q];
00704 }
00705 }
00706 }
00707 }
00708 }
00709 else
00710 {
00711
00712 for (int c=0; c<JVol.numCells(); c++, offset+=nNodes())
00713 {
00714 double detJ = fabs(JVol.detJ()[c]);
00715 int fc = 0;
00716 if (nFacetCases() != 1) fc = facetIndex[c];
00717 const Array<double>& w = W_[fc];
00718 for (int q=0; q<nQuad(); q++, coeffPtr++)
00719 {
00720 double f = (*coeffPtr)*detJ;
00721 for (int n=0; n<nNodes(); n++)
00722 {
00723 aPtr[offset+n] += f*w[n + nNodes()*q];
00724 }
00725 }
00726 }
00727 }
00728 addFlops( JVol.numCells() * (1 + nQuad() * (1 + 2*nNodes())) );
00729 }
00730 else
00731 {
00732 createTwoFormTransformationMatrix(JTrans, JVol);
00733 double* GPtr;
00734
00735 if (testDerivOrder() == 0)
00736 {
00737 GPtr = &(G(beta())[0]);
00738 SUNDANCE_MSG2(transformVerb(),
00739 Tabs() << "transformation matrix=" << G(beta()));
00740 }
00741 else if (unkDerivOrder() == 0)
00742 {
00743 GPtr = &(G(alpha())[0]);
00744 SUNDANCE_MSG2(transformVerb(),
00745 Tabs() << "transformation matrix=" << G(alpha()));
00746 }
00747 else
00748 {
00749 GPtr = &(G(alpha(), beta())[0]);
00750 SUNDANCE_MSG2(transformVerb(),
00751 Tabs() << "transformation matrix="
00752 << G(alpha(),beta()));
00753 }
00754
00755
00756 if (useSumFirstMethod())
00757 {
00758 transformSummingFirst(JTrans.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00759 }
00760 else
00761 {
00762 transformSummingLast(JTrans.numCells(), facetIndex, cellLIDs, GPtr, coeff, A);
00763 }
00764 }
00765 }
00766
00767 void QuadratureIntegral
00768 ::transformSummingFirst(int nCells,
00769 const Array<int>& facetIndex,
00770 const RCP<Array<int> >& cellLIDs,
00771 const double* const GPtr,
00772 const double* const coeff,
00773 RCP<Array<double> >& A) const
00774 {
00775 double* aPtr = &((*A)[0]);
00776 double* coeffPtr = (double*) coeff;
00777 const Array<int>* cellLID = cellLIDs.get();
00778
00779 int transSize = 0;
00780 if (order()==2)
00781 {
00782 transSize = nRefDerivTest() * nRefDerivUnk();
00783 }
00784 else
00785 {
00786 transSize = nRefDerivTest();
00787 }
00788
00789
00790 static Array<double> sumWorkspace;
00791
00792 int swSize = transSize * nNodes();
00793 sumWorkspace.resize(swSize);
00794
00795
00796
00797
00798
00799
00800
00801
00802
00803
00804
00805 if (globalCurve().isCurveValid())
00806 {
00807
00808 Array<double> quadWeightsTmp = quadWeights_[0];
00809 Array<Point> quadPointsTmp = quadPts_[0];
00810 bool isCutByCurve;
00811
00812 for (int c=0; c<nCells; c++)
00813 {
00814
00815 for (int i=0; i<swSize; i++) sumWorkspace[i]=0.0;
00816 int fc = 0;
00817 if (nFacetCases() > 1) fc = facetIndex[c];
00818
00819 const Array<double>& w = W_[fc];
00820
00821 quadWeightsTmp = quadWeights_[fc];
00822 quadPointsTmp = quadPts_[fc];
00823
00824 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00825 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00826 if (isCutByCurve){
00827 Array<double> wi;
00828 if (order()==1){
00829 wi.resize(nQuad() * nRefDerivTest() * nNodesTest());
00830 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00831 for (int n = 0 ; n < nNodes() ; n++)
00832 for (int t=0; t<nRefDerivTest(); t++)
00833 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00834
00835 wi[n + nNodes()*(t + nRefDerivTest()*q)] +=
00836 chop(quadWeightsTmp[q] * W_ACI_F1_[fc][q][t][n]);
00837 }else{
00838 wi.resize(nQuad() * nRefDerivTest() * nNodesTest()
00839 * nRefDerivUnk() * nNodesUnk() );
00840 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00841
00842 for (int t=0; t<nRefDerivTest(); t++)
00843 for (int nt = 0 ; nt < nNodesTest() ; nt++){
00844 for (int u=0; u<nRefDerivUnk(); u++)
00845 for (int nu=0; nu<nNodesUnk(); nu++)
00846 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00847
00848
00849 wi[nu + nNodesUnk()*(nt + nNodesTest()*(u + nRefDerivUnk()*(t + nRefDerivTest()*q)))] +=
00850 chop(quadWeightsTmp[q] * W_ACI_F2_[fc][q][t][nt][u][nu]);
00851 }
00852 }
00853 for (int q=0; q<nQuad(); q++, coeffPtr++){
00854 double f = (*coeffPtr);
00855 for (int n=0; n<swSize; n++){
00856 sumWorkspace[n] += f*wi[n + q*swSize];
00857 }
00858 }
00859 }
00860 else{
00861 for (int q=0; q<nQuad(); q++, coeffPtr++)
00862 {
00863 double f = (*coeffPtr);
00864 for (int n=0; n<swSize; n++)
00865 {
00866 sumWorkspace[n] += f*w[n + q*swSize];
00867 }
00868 }
00869 }
00870 }
00871 }
00872 else
00873 {
00874 for (int c=0; c<nCells; c++)
00875 {
00876
00877 for (int i=0; i<swSize; i++) sumWorkspace[i]=0.0;
00878 int fc = 0;
00879 if (nFacetCases() > 1) fc = facetIndex[c];
00880
00881 const Array<double>& w = W_[fc];
00882
00883 for (int q=0; q<nQuad(); q++, coeffPtr++)
00884 {
00885 double f = (*coeffPtr);
00886 for (int n=0; n<swSize; n++)
00887 {
00888 sumWorkspace[n] += f*w[n + q*swSize];
00889 }
00890 }
00891
00892
00893 const double * const gCell = &(GPtr[transSize*c]);
00894 double* aCell = aPtr + nNodes()*c;
00895 for (int i=0; i<nNodes(); i++)
00896 {
00897 for (int j=0; j<transSize; j++)
00898 {
00899 aCell[i] += sumWorkspace[nNodes()*j + i] * gCell[j];
00900 }
00901 }
00902 }
00903 }
00904 int flops = 2*(nCells * nNodes() * transSize) * (nQuad() + 1) ;
00905 addFlops(flops);
00906 }
00907
00908 void QuadratureIntegral
00909 ::transformSummingLast(int nCells,
00910 const Array<int>& facetIndex,
00911 const RCP<Array<int> >& cellLIDs,
00912 const double* const GPtr,
00913 const double* const coeff,
00914 RCP<Array<double> >& A) const
00915 {
00916 double* aPtr = &((*A)[0]);
00917 int transSize = 0;
00918 const Array<int>* cellLID = cellLIDs.get();
00919
00920 if (order()==2)
00921 {
00922 transSize = nRefDerivTest() * nRefDerivUnk();
00923 }
00924 else
00925 {
00926 transSize = nRefDerivTest();
00927 }
00928
00929
00930
00931 static Array<double> jWorkspace;
00932 jWorkspace.resize(transSize);
00933
00934
00935
00936
00937
00938
00939
00940
00941
00942
00943 if (globalCurve().isCurveValid()){
00944
00945 Array<double> quadWeightsTmp = quadWeights_[0];
00946 Array<Point> quadPointsTmp = quadPts_[0];
00947 bool isCutByCurve;
00948
00949 for (int c=0; c<nCells; c++)
00950 {
00951 const double* const gCell = &(GPtr[transSize*c]);
00952 double* aCell = aPtr + nNodes()*c;
00953 int fc = 0;
00954 if (nFacetCases() > 1) fc = facetIndex[c];
00955 const Array<double>& w = W_[fc];
00956
00957 quadWeightsTmp = quadWeights_[fc];
00958 quadPointsTmp = quadPts_[fc];
00959
00960 quad_.getAdaptedWeights(cellType(), dim(), (*cellLID)[c] , fc ,
00961 mesh() , globalCurve() , quadPointsTmp , quadWeightsTmp , isCutByCurve );
00962 if (isCutByCurve){
00963 Array<double> wi;
00964 if (order()==1){
00965 wi.resize(nQuad() * nRefDerivTest() * nNodesTest());
00966 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00967 for (int n = 0 ; n < nNodes() ; n++)
00968 for (int t=0; t<nRefDerivTest(); t++)
00969 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00970
00971 wi[n + nNodes()*(t + nRefDerivTest()*q)] +=
00972 chop(quadWeightsTmp[q] * W_ACI_F1_[fc][q][t][n]);
00973 }else{
00974 wi.resize(nQuad() * nRefDerivTest() * nNodesTest()
00975 * nRefDerivUnk() * nNodesUnk() );
00976 for (int ii = 0 ; ii < wi.size() ; ii++ ) wi[ii] = 0.0;
00977
00978 for (int t=0; t<nRefDerivTest(); t++)
00979 for (int nt = 0 ; nt < nNodesTest() ; nt++){
00980 for (int u=0; u<nRefDerivUnk(); u++)
00981 for (int nu=0; nu<nNodesUnk(); nu++)
00982 for (int q=0 ; q < quadWeightsTmp.size() ; q++)
00983
00984
00985 wi[nu + nNodesUnk()*(nt + nNodesTest()*(u + nRefDerivUnk()*(t + nRefDerivTest()*q)))] +=
00986 chop(quadWeightsTmp[q] * W_ACI_F2_[fc][q][t][nt][u][nu]);
00987 }
00988 }
00989 for (int q=0; q<nQuad(); q++){
00990 double f = coeff[c*nQuad() + q];
00991 for (int t=0; t<transSize; t++) jWorkspace[t]=f*gCell[t];
00992
00993 for (int n=0; n<nNodes(); n++){
00994 for (int t=0; t<transSize; t++){
00995 aCell[n] += jWorkspace[t]*wi[n + nNodes()*(t + transSize*q)];
00996 }
00997 }
00998 }
00999 }
01000 else{
01001 for (int q=0; q<nQuad(); q++){
01002 double f = coeff[c*nQuad() + q];
01003 for (int t=0; t<transSize; t++) jWorkspace[t]=f*gCell[t];
01004
01005 for (int n=0; n<nNodes(); n++){
01006 for (int t=0; t<transSize; t++){
01007 aCell[n] += jWorkspace[t]*w[n + nNodes()*(t + transSize*q)];
01008 }
01009 }
01010 }
01011 }
01012 }
01013 }
01014 else
01015 {
01016 for (int c=0; c<nCells; c++)
01017 {
01018 const double* const gCell = &(GPtr[transSize*c]);
01019 double* aCell = aPtr + nNodes()*c;
01020 int fc = 0;
01021 if (nFacetCases() > 1) fc = facetIndex[c];
01022 const Array<double>& w = W_[fc];
01023 for (int q=0; q<nQuad(); q++)
01024 {
01025 double f = coeff[c*nQuad() + q];
01026 for (int t=0; t<transSize; t++) jWorkspace[t]=f*gCell[t];
01027
01028 for (int n=0; n<nNodes(); n++)
01029 {
01030 for (int t=0; t<transSize; t++)
01031 {
01032 aCell[n] += jWorkspace[t]*w[n + nNodes()*(t + transSize*q)];
01033 }
01034 }
01035 }
01036 }
01037 }
01038 int flops = nCells * nQuad() * transSize * (1 + 2*nNodes()) ;
01039 addFlops(flops);
01040 }