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