RefIntegral::RefIntegral(int spatialDim, const CellType& maxCellType, int dim, const CellType& cellType, const BasisFamily& testBasis, int alpha, int testDerivOrder, const BasisFamily& unkBasis, int beta, int unkDerivOrder, const QuadratureFamily& quad_in, bool isInternalBdry, const ParametrizedCurve& globalCurve, const Mesh& mesh, int verb) : ElementIntegral(spatialDim, maxCellType, dim, cellType, testBasis, alpha, testDerivOrder, unkBasis, beta, unkDerivOrder, isInternalBdry, globalCurve , mesh ,verb), W_() { Tabs tab0(0); SUNDANCE_MSG1(setupVerb(), tab0 << "************* creating reference 2-form integrals ********"); if (setupVerb()) describe(Out::os()); assertBilinearForm(); W_.resize(nFacetCases()); W_ACI_F2_.resize(nFacetCases()); QuadratureType qType = new GaussianQuadratureType(); int reqOrder = qType.findValidOrder(cellType, std::max(1, unkBasis.order() + testBasis.order())); SUNDANCE_MSG2(setupVerb(), tab0 << "using quadrature order=" << reqOrder); QuadratureFamily quad = qType.createQuadFamily(reqOrder); /* If we have a valid curve (in case of Adaptive Cell Integration) * then we have to choose the quadrature which the user specified*/ if (globalCurve.isCurveValid()){ quad = quad_in; Tabs tab1; SUNDANCE_MSG1(setupVerb(),tab1 << "ACI change quadrature to Quadrature of order: "<<quad.order()); } quad_ = quad; SUNDANCE_MSG2(setupVerb(), tab0 << "processing evaluation cases"); for (int fc=0; fc<nFacetCases(); fc++) { Tabs tab1; SUNDANCE_MSG1(setupVerb(), tab1 << "------ evaluation case " << fc << " of " << nFacetCases() << "-------"); W_[fc].resize(nRefDerivTest() * nNodesTest() * nRefDerivUnk() * nNodesUnk()); for (int i=0; i<W_[fc].size(); i++) W_[fc][i]=0.0; Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest()); Array<Array<Array<Array<double> > > > unkBasisVals(nRefDerivUnk()); getQuad(quad, fc, quadPts_, quadWeights_); int nQuad = quadPts_.size(); for (int r=0; r<nRefDerivTest(); r++) { Tabs tab2; SUNDANCE_MSG2(setupVerb(), tab2 << "evaluating test function basis derivative " << r << " of " << nRefDerivTest()); testBasisVals[r].resize(testBasis.dim()); MultiIndex mi; if (testDerivOrder==1) mi[r] = 1; SpatialDerivSpecifier deriv(mi); testBasis.refEval(evalCellType(), quadPts_, deriv, testBasisVals[r], setupVerb()); } for (int r=0; r<nRefDerivUnk(); r++) { Tabs tab2; SUNDANCE_MSG2(setupVerb(), tab2 << "evaluating unknown function basis derivative " << r << " of " << nRefDerivUnk()); unkBasisVals[r].resize(unkBasis.dim()); MultiIndex mi; if (unkDerivOrder==1) mi[r] = 1; SpatialDerivSpecifier deriv(mi); unkBasis.refEval(evalCellType(), quadPts_, deriv, unkBasisVals[r], setupVerb()); } SUNDANCE_MSG2(setupVerb(), tab1 << "doing quadrature..."); int vecComp = 0; W_ACI_F2_[fc].resize(nQuad); for (int q=0; q<nQuad; q++) { W_ACI_F2_[fc][q].resize(nRefDerivTest()); for (int t=0; t<nRefDerivTest(); t++) { W_ACI_F2_[fc][q][t].resize(nNodesTest()); for (int nt=0; nt<nNodesTest(); nt++) { W_ACI_F2_[fc][q][t][nt].resize(nRefDerivUnk()); for (int u=0; u<nRefDerivUnk(); u++) { W_ACI_F2_[fc][q][t][nt][u].resize(nNodesUnk()); for (int nu=0; nu<nNodesUnk(); nu++) { value(fc, t, nt, u, nu) += chop(quadWeights_[q] * testBasisVals[t][vecComp][q][nt] * unkBasisVals[u][vecComp][q][nu]); W_ACI_F2_[fc][q][t][nt][u][nu] = chop( testBasisVals[t][vecComp][q][nt] * unkBasisVals[u][vecComp][q][nu] ); } } } } } SUNDANCE_MSG2(setupVerb(), tab1 << "...done"); addFlops(4*nQuad*nRefDerivTest()*nNodesTest()*nRefDerivUnk()*nNodesUnk() + W_[fc].size()); for (int i=0; i<W_[fc].size(); i++) W_[fc][i] = chop(W_[fc][i]); } SUNDANCE_MSG1(setupVerb(), tab0 << "----------------------------------------"); SUNDANCE_MSG4(setupVerb(), tab0 << "reference bilinear form integral results"); if (setupVerb() >= 4) { for (int fc=0; fc<nFacetCases(); fc++) { Tabs tab1; SUNDANCE_MSG4(setupVerb(), tab1 << "evaluation case " << fc << " of " << nFacetCases()); for (int rt=0; rt<nRefDerivTest(); rt++) { for (int ru=0; ru<nRefDerivUnk(); ru++) { Tabs tab2; MultiIndex miTest; if (testDerivOrder==1) miTest[rt] = 1; MultiIndex miUnk; if (unkDerivOrder==1) miUnk[ru] = 1; SUNDANCE_MSG1(setupVerb(), tab2 << "test multiindex=" << miTest << " unk multiindex=" << miUnk); ios_base::fmtflags oldFlags = Out::os().flags(); Out::os().setf(ios_base::right); Out::os().setf(ios_base::showpoint); for (int nt=0; nt<nNodesTest(); nt++) { Tabs tab3; Out::os() << tab3 << setw(10) << nt; for (int nu=0; nu<nNodesUnk(); nu++) { Out::os() << setw(12) << std::setprecision(5) << value(fc, rt, nt, ru, nu) ; } Out::os() << std::endl; } Out::os().flags(oldFlags); } } } } SUNDANCE_MSG1(setupVerb(), tab0 << "done reference bilinear form ctor"); }
RefIntegral::RefIntegral(int spatialDim, const CellType& maxCellType, int dim, const CellType& cellType, const QuadratureFamily& quad_in, bool isInternalBdry, const ParametrizedCurve& globalCurve, const Mesh& mesh, int verb) : ElementIntegral(spatialDim, maxCellType, dim, cellType, isInternalBdry, globalCurve , mesh, verb), W_() { Tabs tab0(0); SUNDANCE_MSG1(setupVerb(), tab0 << "************* creating reference 0-form integrals ********"); if (setupVerb()) describe(Out::os()); /* we need to sum the quadrature weights to compute the volume of the reference cell */ QuadratureFamily quad = new GaussianQuadrature(2); /* If we have a valid curve (in case of Adaptive Cell Integration) * then we have to choose the quadrature which the user specified*/ if (globalCurve.isCurveValid()){ quad = quad_in; Tabs tab1; SUNDANCE_MSG1(setupVerb(),tab1 << "ACI change quadrature to Quadrature of order: "<<quad.order()); } quad_ = quad; Array<Point> quadPts; Array<double> quadWeights; W_.resize(1); W_[0].resize(1); quad.getPoints(cellType, quadPts, quadWeights); quadWeights_ = quadWeights; for (int q=0; q<quadWeights.size(); q++) { W_[0][0] += quadWeights[q]; } }
RefIntegral::RefIntegral(int spatialDim, const CellType& maxCellType, int dim, const CellType& cellType, const BasisFamily& testBasis, int alpha, int testDerivOrder, const QuadratureFamily& quad_in, bool isInternalBdry, const ParametrizedCurve& globalCurve, const Mesh& mesh, int verb) : ElementIntegral(spatialDim, maxCellType, dim, cellType, testBasis, alpha, testDerivOrder, isInternalBdry, globalCurve , mesh , verb), W_() { Tabs tab0(0); SUNDANCE_MSG1(setupVerb(), tab0 << "************* creating reference 1-form integrals ********"); if (setupVerb()) describe(Out::os()); assertLinearForm(); W_.resize(nFacetCases()); W_ACI_F1_.resize(nFacetCases()); /* Determine the quadrature order needed for exact integrations */ QuadratureType qType = new GaussianQuadratureType(); int reqOrder = qType.findValidOrder(cellType, std::max(1, testBasis.order())); SUNDANCE_MSG2(setupVerb(), tab0 << "using quadrature order=" << reqOrder); /* Create a quadrature family of the required order */ QuadratureFamily quad = qType.createQuadFamily(reqOrder); /* If we have a valid curve (in case of Adaptive Cell Integration) * then we have to choose the quadrature which the user specified*/ if (globalCurve.isCurveValid()){ quad = quad_in; Tabs tab1; SUNDANCE_MSG1(setupVerb(),tab1 << "ACI change quadrature to Quadrature of order: "<<quad.order()); } quad_ = quad; /* We now loop over the different evaluation cases, integrating the * basis functions for each. Because this is a reference integral, * we can actually do the untransformed integrals here. */ for (int fc=0; fc<nFacetCases(); fc++) { Tabs tab1; SUNDANCE_MSG2(setupVerb(), tab1 << "evaluation case=" << fc << " of " << nFacetCases()); /* initialize size of untransformed integral results array */ W_[fc].resize(nRefDerivTest() * nNodesTest()); /* initialize values of integrals to zero */ for (int i=0; i<W_[fc].size(); i++) { W_[fc][i]=0.0; } Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest()); /* get quadrature points */ getQuad(quad, fc, quadPts_, quadWeights_); int nQuad = quadPts_.size(); /* compute the basis functions */ for (int r=0; r<nRefDerivTest(); r++) { Tabs tab2; SUNDANCE_MSG2(setupVerb(), tab2 << "evaluating basis derivative " << r << " of " << nRefDerivTest()); testBasisVals[r].resize(testBasis.dim()); MultiIndex mi; if (testDerivOrder==1) mi[r] = 1; SpatialDerivSpecifier deriv(mi); testBasis.refEval(evalCellType(), quadPts_, deriv, testBasisVals[r], setupVerb()); } /* do the quadrature */ SUNDANCE_MSG2(setupVerb(), tab1 << "doing quadrature"); int vecComp = 0; W_ACI_F1_[fc].resize(nQuad); for (int q=0; q<nQuad; q++) { W_ACI_F1_[fc][q].resize(nRefDerivTest()); for (int t=0; t<nRefDerivTest(); t++) { W_ACI_F1_[fc][q][t].resize(nNodesTest()); for (int nt=0; nt<nNodesTest(); nt++) { value(fc, t, nt) += chop(quadWeights_[q] * testBasisVals[t][vecComp][q][nt]) ; W_ACI_F1_[fc][q][t][nt] = chop(testBasisVals[t][vecComp][q][nt]); } } } for (int i=0; i<W_[fc].size(); i++) W_[fc][i] = chop(W_[fc][i]); addFlops(3*nQuad*nRefDerivTest()*nNodesTest() + W_[fc].size()); } /* print the result */ SUNDANCE_MSG4(setupVerb(), tab0 << "--------------------------------------"); SUNDANCE_MSG4(setupVerb(), tab0 << "reference linear form integral results"); if (setupVerb() >= 4) { for (int fc=0; fc<nFacetCases(); fc++) { Tabs tab1; SUNDANCE_MSG4(setupVerb(), tab1 << "------ evaluation case " << fc << " of " << nFacetCases() << "-------"); for (int r=0; r<nRefDerivTest(); r++) { Tabs tab2; MultiIndex mi; if (testDerivOrder==1) mi[r] = 1; SUNDANCE_MSG1(setupVerb(), tab2 << "multiindex=" << mi); ios_base::fmtflags oldFlags = Out::os().flags(); Out::os().setf(ios_base::right); Out::os().setf(ios_base::showpoint); for (int nt=0; nt<nNodesTest(); nt++) { Tabs tab3; Out::os() << tab3 << setw(10) << nt << setw(12) << std::setprecision(5) << value(fc, r, nt) << std::endl; } Out::os().flags(oldFlags); } } } SUNDANCE_MSG1(setupVerb(), tab0 << "done reference linear form ctor"); }