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 "SundanceReducedIntegral.hpp"
00043 #include "SundanceGaussianQuadrature.hpp"
00044 #include "SundanceGaussianQuadratureType.hpp"
00045 #include "SundanceQuadratureType.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::ios_base;
00055 using std::setw;
00056 using std::endl;
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& reduced0IntegrationTimer()
00070 {
00071 static RCP<Time> rtn
00072 = TimeMonitor::getNewTimer("reduced 0-form integration");
00073 return *rtn;
00074 }
00075
00076
00077 static Time& reduced1IntegrationTimer()
00078 {
00079 static RCP<Time> rtn
00080 = TimeMonitor::getNewTimer("reduced 1-form integration");
00081 return *rtn;
00082 }
00083
00084
00085 static Time& reduced2IntegrationTimer()
00086 {
00087 static RCP<Time> rtn
00088 = TimeMonitor::getNewTimer("reduced 2-form integration");
00089 return *rtn;
00090 }
00091
00092
00093
00094 ReducedIntegral::ReducedIntegral(int spatialDim,
00095 const CellType& maxCellType,
00096 int dim,
00097 const CellType& cellType,
00098 const QuadratureFamily& quad_in,
00099 bool isInternalBdry,
00100 const ParametrizedCurve& globalCurve,
00101 const Mesh& mesh,
00102 int verb)
00103 : ElementIntegral(spatialDim, maxCellType, dim, cellType, isInternalBdry, globalCurve , mesh,
00104 verb), W_()
00105 {
00106 Tabs tab0(0);
00107
00108 SUNDANCE_MSG1(setupVerb(),
00109 tab0 << "************* creating reduced 0-form integrals ********");
00110 if (setupVerb()) describe(Out::os());
00111
00112
00113
00114 QuadratureFamily quad = new GaussianQuadrature(1);
00115
00116
00117 TEUCHOS_TEST_FOR_EXCEPT(globalCurve.isCurveValid());
00118
00119 Array<Point> quadPts;
00120 Array<double> quadWeights;
00121
00122 W_.resize(1);
00123 W_[0].resize(1);
00124
00125 quad.getPoints(cellType, quadPts, quadWeights);
00126
00127 for (int q=0; q<quadWeights.size(); q++) {
00128 W_[0][0] += quadWeights[q];
00129 }
00130 }
00131
00132 ReducedIntegral::ReducedIntegral(int spatialDim,
00133 const CellType& maxCellType,
00134 int dim,
00135 const CellType& cellType,
00136 const BasisFamily& testBasis,
00137 int alpha,
00138 int testDerivOrder,
00139 const QuadratureFamily& quad_in,
00140 bool isInternalBdry,
00141 const ParametrizedCurve& globalCurve,
00142 const Mesh& mesh,
00143 int verb)
00144 : ElementIntegral(spatialDim, maxCellType, dim, cellType,
00145 testBasis, alpha, testDerivOrder, isInternalBdry, globalCurve , mesh , verb), W_()
00146 {
00147 Tabs tab0(0);
00148 SUNDANCE_MSG1(setupVerb(),
00149 tab0 << "************* creating reduced 1-form integrals ********");
00150 if (setupVerb()) describe(Out::os());
00151 assertLinearForm();
00152
00153 W_.resize(nFacetCases());
00154
00155
00156 QuadratureType qType = new GaussianQuadratureType();
00157 int reqOrder = qType.findValidOrder(cellType,
00158 std::max(1, testBasis.order()));
00159
00160 SUNDANCE_MSG2(setupVerb(),
00161 tab0 << "using quadrature order=" << reqOrder);
00162
00163
00164 QuadratureFamily quad = qType.createQuadFamily(reqOrder);
00165
00166
00167 TEUCHOS_TEST_FOR_EXCEPT(globalCurve.isCurveValid());
00168
00169
00170
00171
00172 for (int fc=0; fc<nFacetCases(); fc++)
00173 {
00174 Tabs tab1;
00175 SUNDANCE_MSG2(setupVerb(),
00176 tab1 << "evaluation case=" << fc << " of " << nFacetCases());
00177
00178 W_[fc].resize(nRefDerivTest() * nNodesTest());
00179
00180
00181 for (int i=0; i<W_[fc].size(); i++) { W_[fc][i]=0.0; }
00182
00183 Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());
00184
00185
00186
00187 Array<Point> quadPts;
00188 Array<double> quadWeights;
00189 getQuad(quad, fc, quadPts, quadWeights);
00190
00191 int nQuad = quadPts.size();
00192
00193
00194 for (int r=0; r<nRefDerivTest(); r++)
00195 {
00196 Tabs tab2;
00197 SUNDANCE_MSG2(setupVerb(), tab2 << "evaluating basis derivative "
00198 << r << " of " << nRefDerivTest());
00199
00200 testBasisVals[r].resize(testBasis.dim());
00201 MultiIndex mi;
00202 if (testDerivOrder==1) mi[r] = 1;
00203 SpatialDerivSpecifier deriv(mi);
00204 testBasis.refEval(evalCellType(), quadPts, deriv,
00205 testBasisVals[r], setupVerb());
00206 }
00207
00208
00209 SUNDANCE_MSG2(setupVerb(), tab1 << "doing quadrature");
00210 int vecComp = 0;
00211 for (int q=0; q<nQuad; q++)
00212 {
00213 for (int t=0; t<nRefDerivTest(); t++)
00214 {
00215 for (int nt=0; nt<nNodesTest(); nt++)
00216 {
00217 value(fc, t, nt)
00218 += chop(quadWeights[q] * testBasisVals[t][vecComp][q][nt]) ;
00219 }
00220 }
00221 }
00222
00223 for (int i=0; i<W_[fc].size(); i++) W_[fc][i] = chop(W_[fc][i]);
00224
00225 addFlops(3*nQuad*nRefDerivTest()*nNodesTest() + W_[fc].size());
00226 }
00227
00228
00229 SUNDANCE_MSG4(setupVerb(), tab0 << "--------------------------------------");
00230 SUNDANCE_MSG4(setupVerb(), tab0 << "reduced linear form integral results");
00231 if (setupVerb() >= 4)
00232 {
00233 for (int fc=0; fc<nFacetCases(); fc++)
00234 {
00235 Tabs tab1;
00236 SUNDANCE_MSG4(setupVerb(), tab1 << "------ evaluation case " << fc << " of "
00237 << nFacetCases() << "-------");
00238
00239 for (int r=0; r<nRefDerivTest(); r++)
00240 {
00241 Tabs tab2;
00242
00243 MultiIndex mi;
00244 if (testDerivOrder==1) mi[r] = 1;
00245 SUNDANCE_MSG1(setupVerb(), tab2 << "multiindex=" << mi);
00246
00247 ios_base::fmtflags oldFlags = Out::os().flags();
00248 Out::os().setf(ios_base::right);
00249 Out::os().setf(ios_base::showpoint);
00250 for (int nt=0; nt<nNodesTest(); nt++)
00251 {
00252 Tabs tab3;
00253 Out::os() << tab3 << setw(10) << nt
00254 << setw(12) << std::setprecision(5) << value(fc, r, nt)
00255 << std::endl;
00256 }
00257 Out::os().flags(oldFlags);
00258 }
00259 }
00260 }
00261
00262 SUNDANCE_MSG1(setupVerb(), tab0 << "done reduced linear form ctor");
00263 }
00264
00265
00266
00267
00268 ReducedIntegral::ReducedIntegral(int spatialDim,
00269 const CellType& maxCellType,
00270 int dim,
00271 const CellType& cellType,
00272 const BasisFamily& testBasis,
00273 int alpha,
00274 int testDerivOrder,
00275 const BasisFamily& unkBasis,
00276 int beta,
00277 int unkDerivOrder,
00278 const QuadratureFamily& quad_in,
00279 bool isInternalBdry,
00280 const ParametrizedCurve& globalCurve,
00281 const Mesh& mesh,
00282 int verb)
00283 : ElementIntegral(spatialDim, maxCellType, dim, cellType,
00284 testBasis, alpha, testDerivOrder,
00285 unkBasis, beta, unkDerivOrder, isInternalBdry, globalCurve , mesh ,verb), W_()
00286
00287 {
00288 Tabs tab0(0);
00289 SUNDANCE_MSG1(setupVerb(),
00290 tab0 << "************* creating reduced 2-form integrals ********");
00291 if (setupVerb()) describe(Out::os());
00292
00293 assertBilinearForm();
00294
00295 W_.resize(nFacetCases());
00296
00297 QuadratureType qType = new GaussianQuadratureType();
00298 int reqOrder = qType.findValidOrder(cellType,
00299 std::max(1, unkBasis.order() + testBasis.order()));
00300
00301 SUNDANCE_MSG2(setupVerb(),
00302 tab0 << "using quadrature order=" << reqOrder);
00303 QuadratureFamily quad = qType.createQuadFamily(reqOrder);
00304
00305
00306 TEUCHOS_TEST_FOR_EXCEPT(globalCurve.isCurveValid());
00307
00308
00309 SUNDANCE_MSG2(setupVerb(),
00310 tab0 << "processing evaluation cases");
00311
00312 for (int fc=0; fc<nFacetCases(); fc++)
00313 {
00314 Tabs tab1;
00315 SUNDANCE_MSG1(setupVerb(), tab1 << "------ evaluation case " << fc << " of "
00316 << nFacetCases() << "-------");
00317
00318 W_[fc].resize(nRefDerivTest() * nNodesTest() * nRefDerivUnk() * nNodesUnk());
00319 for (int i=0; i<W_[fc].size(); i++) W_[fc][i]=0.0;
00320
00321 Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());
00322 Array<Array<Array<Array<double> > > > unkBasisVals(nRefDerivUnk());
00323
00324 Array<Point> quadPts;
00325 Array<double> quadWeights;
00326 getQuad(quad, fc, quadPts, quadWeights);
00327
00328 int nQuad = quadPts.size();
00329
00330 for (int r=0; r<nRefDerivTest(); r++)
00331 {
00332 Tabs tab2;
00333 SUNDANCE_MSG2(setupVerb(), tab2
00334 << "evaluating test function basis derivative "
00335 << r << " of " << nRefDerivTest());
00336 testBasisVals[r].resize(testBasis.dim());
00337 MultiIndex mi;
00338 if (testDerivOrder==1) mi[r] = 1;
00339 SpatialDerivSpecifier deriv(mi);
00340 testBasis.refEval(evalCellType(), quadPts, deriv,
00341 testBasisVals[r], setupVerb());
00342 }
00343
00344 for (int r=0; r<nRefDerivUnk(); r++)
00345 {
00346 Tabs tab2;
00347 SUNDANCE_MSG2(setupVerb(), tab2
00348 << "evaluating unknown function basis derivative "
00349 << r << " of " << nRefDerivUnk());
00350 unkBasisVals[r].resize(unkBasis.dim());
00351 MultiIndex mi;
00352 if (unkDerivOrder==1) mi[r] = 1;
00353 SpatialDerivSpecifier deriv(mi);
00354 unkBasis.refEval(evalCellType(),
00355 quadPts, deriv, unkBasisVals[r], setupVerb());
00356 }
00357
00358 SUNDANCE_MSG2(setupVerb(), tab1 << "doing quadrature...");
00359 int vecComp = 0;
00360 for (int q=0; q<nQuad; q++)
00361 {
00362 for (int t=0; t<nRefDerivTest(); t++)
00363 {
00364 for (int nt=0; nt<nNodesTest(); nt++)
00365 {
00366 for (int u=0; u<nRefDerivUnk(); u++)
00367 {
00368 for (int nu=0; nu<nNodesUnk(); nu++)
00369 {
00370 value(fc, t, nt, u, nu)
00371 += chop(quadWeights[q] * testBasisVals[t][vecComp][q][nt]
00372 * unkBasisVals[u][vecComp][q][nu]);
00373 }
00374 }
00375 }
00376 }
00377 }
00378 SUNDANCE_MSG2(setupVerb(), tab1 << "...done");
00379 addFlops(4*nQuad*nRefDerivTest()*nNodesTest()*nRefDerivUnk()*nNodesUnk()
00380 + W_[fc].size());
00381 for (int i=0; i<W_[fc].size(); i++) W_[fc][i] = chop(W_[fc][i]);
00382 }
00383
00384 SUNDANCE_MSG1(setupVerb(), tab0
00385 << "----------------------------------------");
00386 SUNDANCE_MSG4(setupVerb(), tab0
00387 << "reduced bilinear form integral results");
00388 if (setupVerb() >= 4)
00389 {
00390 for (int fc=0; fc<nFacetCases(); fc++)
00391 {
00392 Tabs tab1;
00393 SUNDANCE_MSG4(setupVerb(), tab1 << "evaluation case " << fc << " of "
00394 << nFacetCases());
00395
00396 for (int rt=0; rt<nRefDerivTest(); rt++)
00397 {
00398 for (int ru=0; ru<nRefDerivUnk(); ru++)
00399 {
00400 Tabs tab2;
00401 MultiIndex miTest;
00402 if (testDerivOrder==1) miTest[rt] = 1;
00403 MultiIndex miUnk;
00404 if (unkDerivOrder==1) miUnk[ru] = 1;
00405 SUNDANCE_MSG1(setupVerb(), tab2 << "test multiindex=" << miTest
00406 << " unk multiindex=" << miUnk);
00407
00408 ios_base::fmtflags oldFlags = Out::os().flags();
00409 Out::os().setf(ios_base::right);
00410 Out::os().setf(ios_base::showpoint);
00411 for (int nt=0; nt<nNodesTest(); nt++)
00412 {
00413 Tabs tab3;
00414 Out::os() << tab3 << setw(10) << nt;
00415 for (int nu=0; nu<nNodesUnk(); nu++)
00416 {
00417 Out::os() << setw(12) << std::setprecision(5)
00418 << value(fc, rt, nt, ru, nu) ;
00419 }
00420 Out::os() << std::endl;
00421 }
00422 Out::os().flags(oldFlags);
00423 }
00424 }
00425 }
00426 }
00427
00428 SUNDANCE_MSG1(setupVerb(), tab0 << "done reduced bilinear form ctor");
00429 }
00430
00431
00432
00433
00434 void ReducedIntegral::transformZeroForm(const CellJacobianBatch& JTrans,
00435 const CellJacobianBatch& JVol,
00436 const Array<int>& isLocalFlag,
00437 const Array<int>& facetIndex,
00438 const RCP<Array<int> >& cellLIDs,
00439 const double* const coeffs,
00440 RCP<Array<double> >& A) const
00441 {
00442 TimeMonitor timer(reduced0IntegrationTimer());
00443
00444 TEUCHOS_TEST_FOR_EXCEPTION(order() != 0, std::logic_error,
00445 "ReducedIntegral::transformZeroForm() called "
00446 "for form of order " << order());
00447
00448 Tabs tabs;
00449 SUNDANCE_MSG1(integrationVerb(), tabs << "doing zero form by reduced integration");
00450
00451 double& a = (*A)[0];
00452 int flops = 0;
00453
00454
00455
00456
00457 double w = W_[0][0];
00458 if ((int) isLocalFlag.size()==0)
00459 {
00460
00461 if (globalCurve().isCurveValid())
00462 {
00463 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00464 }
00465 else
00466 {
00467 for (int c=0; c<JVol.numCells(); c++)
00468 {
00469 flops+=3;
00470 a += w * coeffs[c] * fabs(JVol.detJ()[c]);
00471 }
00472 }
00473
00474 }
00475 else
00476 {
00477 TEUCHOS_TEST_FOR_EXCEPTION( (int) isLocalFlag.size() != JVol.numCells(),
00478 std::runtime_error,
00479 "mismatch between isLocalFlag.size()="
00480 << isLocalFlag.size()
00481 << " and JVol.numCells()=" << JVol.numCells());
00482
00483 if (globalCurve().isCurveValid())
00484 {
00485 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00486 }
00487 else
00488 {
00489 for (int c=0; c<JVol.numCells(); c++)
00490 {
00491 if (isLocalFlag[c])
00492 {
00493 flops+=3;
00494 a += w * coeffs[c] * fabs(JVol.detJ()[c]);
00495 }
00496 }
00497 }
00498 }
00499 addFlops(flops);
00500 }
00501
00502 void ReducedIntegral::transformOneForm(const CellJacobianBatch& JTrans,
00503 const CellJacobianBatch& JVol,
00504 const Array<int>& facetIndex,
00505 const RCP<Array<int> >& cellLIDs,
00506 const double* const coeffs,
00507 RCP<Array<double> >& A) const
00508 {
00509 TimeMonitor timer(reduced1IntegrationTimer());
00510 TEUCHOS_TEST_FOR_EXCEPTION(order() != 1, std::logic_error,
00511 "ReducedIntegral::transformOneForm() called for form "
00512 "of order " << order());
00513 Tabs tabs;
00514 SUNDANCE_MSG1(integrationVerb(),
00515 tabs << "doing one form by reduced integration");
00516
00517 int nfc = nFacetCases();
00518
00519
00520
00521
00522 if (testDerivOrder() == 0)
00523 {
00524 double* aPtr = &((*A)[0]);
00525 int count = 0;
00526 if (globalCurve().isCurveValid())
00527 {
00528 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00529 }
00530 else
00531 {
00532 for (int c=0; c<JVol.numCells(); c++)
00533 {
00534 double w_cell = coeffs[c] * fabs(JVol.detJ()[c]);
00535 int fc = 0;
00536 if (nfc != 1) fc = facetIndex[c];
00537
00538 const Array<double>& w = W_[fc];
00539 for (int n=0; n<nNodes(); n++, count++)
00540 {
00541 aPtr[count] += w_cell*w[n];
00542 }
00543 }
00544 }
00545 addFlops(JVol.numCells() * (nNodes() + 2));
00546 }
00547 else
00548 {
00549
00550
00551
00552 int nCells = JVol.numCells();
00553 double one = 1.0;
00554 int nTransRows = nRefDerivTest();
00555
00556 createOneFormTransformationMatrix(JTrans, JVol);
00557
00558 SUNDANCE_MSG3(transformVerb(),
00559 Tabs() << "transformation matrix=" << G(alpha()));
00560 int nNodes0 = nNodes();
00561
00562 if (nFacetCases()==1)
00563 {
00564
00565
00566
00567 if (globalCurve().isCurveValid())
00568 {
00569 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00570 }
00571 else
00572 {
00573 int N = 1;
00574 for (int c=0; c<JVol.numCells(); c++)
00575 {
00576 double* aPtr = &((*A)[c*nNodes0]);
00577 ::dgemm_("N", "N", &nNodes0, &N, &nTransRows, &(coeffs[c]), &(W_[0][0]),
00578 &nNodes0, &(G(alpha())[c*nTransRows]), &nTransRows, &one,
00579 aPtr, &nNodes0);
00580 }
00581 }
00582 }
00583 else
00584 {
00585
00586
00587 if (globalCurve().isCurveValid())
00588 {
00589 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00590 }
00591 else
00592 {
00593 int N = 1;
00594 for (int c=0; c<JVol.numCells(); c++)
00595 {
00596 int fc = 0;
00597 if (nfc != 1) fc = facetIndex[c];
00598 double* aPtr = &((*A)[c*nNodes0]);
00599 ::dgemm_("N", "N", &nNodes0, &N, &nTransRows, &(coeffs[c]), &(W_[fc][0]),
00600 &nNodes0, &(G(alpha())[c*nTransRows]), &nTransRows, &one,
00601 aPtr, &nNodes0);
00602 }
00603 }
00604 }
00605
00606 addFlops(2 * nNodes0 * nCells * nTransRows);
00607 }
00608 }
00609
00610 void ReducedIntegral::transformTwoForm(const CellJacobianBatch& JTrans,
00611 const CellJacobianBatch& JVol,
00612 const Array<int>& facetIndex,
00613 const RCP<Array<int> >& cellLIDs,
00614 const double* const coeffs,
00615 RCP<Array<double> >& A) const
00616 {
00617 TimeMonitor timer(reduced2IntegrationTimer());
00618 TEUCHOS_TEST_FOR_EXCEPTION(order() != 2, std::logic_error,
00619 "ReducedIntegral::transformTwoForm() called for form "
00620 "of order " << order());
00621
00622 Tabs tabs;
00623 SUNDANCE_MSG1(transformVerb(), tabs << "doing two form by reduced integration");
00624
00625 int nfc = nFacetCases();
00626
00627
00628
00629 if (testDerivOrder() == 0 && unkDerivOrder() == 0)
00630 {
00631 if (globalCurve().isCurveValid())
00632 {
00633 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00634 }
00635 else
00636 {
00637 double* aPtr = &((*A)[0]);
00638 int count = 0;
00639 for (int c=0; c<JVol.numCells(); c++)
00640 {
00641 int fc = 0;
00642 if (nFacetCases() != 1) fc = facetIndex[c];
00643
00644 const Array<double>& w = W_[fc];
00645 double w_cell = coeffs[c] * fabs(JVol.detJ()[c]);
00646 for (int n=0; n<nNodes(); n++, count++)
00647 {
00648 aPtr[count] += w_cell*w[n];
00649 }
00650 }
00651 }
00652 addFlops(JVol.numCells() * (nNodes() + 1));
00653 }
00654 else
00655 {
00656
00657
00658
00659 int nCells = JVol.numCells();
00660 double one = 1.0;
00661 int nTransRows = nRefDerivUnk()*nRefDerivTest();
00662
00663 createTwoFormTransformationMatrix(JTrans, JVol);
00664
00665 double* GPtr;
00666 if (testDerivOrder() == 0)
00667 {
00668 GPtr = &(G(beta())[0]);
00669 SUNDANCE_MSG2(transformVerb(),
00670 Tabs() << "transformation matrix=" << G(beta()));
00671 }
00672 else if (unkDerivOrder() == 0)
00673 {
00674 GPtr = &(G(alpha())[0]);
00675 SUNDANCE_MSG2(transformVerb(),
00676 Tabs() << "transformation matrix=" << G(alpha()));
00677 }
00678 else
00679 {
00680 GPtr = &(G(alpha(), beta())[0]);
00681 SUNDANCE_MSG2(transformVerb(),
00682 Tabs() << "transformation matrix="
00683 << G(alpha(),beta()));
00684 }
00685
00686 int nNodes0 = nNodes();
00687
00688 if (nFacetCases()==1)
00689 {
00690
00691
00692
00693 if (globalCurve().isCurveValid())
00694 {
00695 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00696 }
00697 else
00698 {
00699 int N = 1;
00700 for (int c=0; c<JVol.numCells(); c++)
00701 {
00702 double* aPtr = &((*A)[c*nNodes0]);
00703 double* gPtr = &(GPtr[c*nTransRows]);
00704 SUNDANCE_MSG2(integrationVerb(),
00705 tabs << "transforming c=" << c << ", W=" << W_[0]);
00706
00707 ::dgemm_("N", "N", &nNodes0, &N, &nTransRows, &(coeffs[c]), &(W_[0][0]),
00708 &nNodes0, gPtr, &nTransRows, &one,
00709 aPtr, &nNodes0);
00710 }
00711 }
00712 }
00713 else
00714 {
00715
00716
00717 if (globalCurve().isCurveValid())
00718 {
00719 TEUCHOS_TEST_FOR_EXCEPT(globalCurve().isCurveValid());
00720 }
00721 else
00722 {
00723 int N = 1;
00724 for (int c=0; c<JVol.numCells(); c++)
00725 {
00726 int fc = 0;
00727 if (nfc != 1) fc = facetIndex[c];
00728 double* aPtr = &((*A)[c*nNodes0]);
00729 double* gPtr = &(GPtr[c*nTransRows]);
00730 SUNDANCE_MSG2(integrationVerb(),
00731 tabs << "c=" << c << ", facet case=" << fc
00732 << " W=" << W_[fc]);
00733
00734 ::dgemm_("N", "N", &nNodes0, &N, &nTransRows, &(coeffs[c]), &(W_[fc][0]),
00735 &nNodes0, gPtr, &nTransRows, &one,
00736 aPtr, &nNodes0);
00737 }
00738 }
00739 }
00740
00741 addFlops(2 * nNodes0 * nCells * nTransRows);
00742 }
00743 }
00744