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