Ejemplo n.º 1
0
void ElementIntegral::getQuad(const QuadratureFamily& quad,
                              int evalCase, Array<Point>& quadPts, Array<double>& quadWeights) const
{
    Tabs tab(0);

    SUNDANCE_MSG2(setupVerb(), tab << "getting quad points for rule "
                  << quad);

    if (nFacetCases()==1)
    {
        quad.getPoints(cellType(), quadPts, quadWeights);
    }
    else
    {
        int dim = dimension(cellType());
        quad.getFacetPoints(maxCellType(), dim, evalCase,
                            quadPts, quadWeights);
    }

    if (setupVerb() >= 4)
    {
        Tabs tab1;
        SUNDANCE_MSG4(setupVerb(),
                      tab1 << "quadrature points on ref cell are:");
        printQuad(Out::os(), quadPts, quadWeights);
    }
}
IQI_HdivLF_DivV_Cell::IQI_HdivLF_DivV_Cell( int spatialDim,
					    const CellType& maxCellType,
					    const BasisFamily& testBasis,
					    const QuadratureFamily& quad,
					    const ParameterList& verbParams):
  ElementIntegralLinearFormCell( spatialDim ,
				 maxCellType ,
				 testBasis ,
				 quad ,
				 verbParams ),
  DivV_( testBasis.nReferenceDOFs(maxCellType,maxCellType) , quad.getNumPoints(maxCellType) ),
  QP_( quad.getNumPoints( maxCellType ) , spatialDim ),
  QW_( quad.getNumPoints( maxCellType ) )
{
  // bypass testBasis.refEval and use Intrepid to fill DivV
  // warning: only works on tets right now.
  Intrepid::Basis_HDIV_TET_I1_FEM<double,Intrepid::FieldContainer<double> > myBasis;

  // move quadrature points into a field container
  Array<Point> qpSundance;
  Array<double> qwSundance;
  quad.getPoints( maxCellType , qpSundance , qwSundance );

  for (int i=0;i<(int)qpSundance.size();i++) {
    for (int j=0;j<3;j++) {
      QP_(i,j) = qpSundance[i][j];
    }
    QW_(i)=qwSundance[i];
  }

  // now tabulate the divergences 
  myBasis.getValues( DivV_ , QP_ , Intrepid::OPERATOR_DIV );

}
MaximalQuadratureIntegral::MaximalQuadratureIntegral(
  const CellType& cellType,
  const QuadratureFamily& quad,
  const ParametrizedCurve& globalCurve,
  const Mesh& mesh,
  int verb)
  : ElementIntegral(dimension(cellType), cellType, dimension(cellType),
    cellType, true, globalCurve, mesh, verb),
    quad_(quad),
    quadPts_(),
    quadWeights_(),
    W_(),
    useSumFirstMethod_(true)
{
  Tabs tab0(0);
  
  SUNDANCE_MSG1(setupVerb(), tab0 << "MaximalQuadratureIntegral ctor for 0-form");
  if (setupVerb()) describe(Out::os());
  
  SUNDANCE_MSG1(setupVerb(), tab0 << "quadrature family=" << quad);  

  /* create the quad points and weights */
  quad.getPoints(cellType, quadPts_, quadWeights_);
  int nQuad = quadPts_.size();
      
  W_.resize(nQuad);
      
  for (int q=0; q<nQuad; q++)
  {
    W_[q] = quadWeights_[q];
  }
}
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];
  }
}
Ejemplo n.º 5
0
CurveEvalMediator::CurveEvalMediator(
		const Mesh& mesh, 	const ParametrizedCurve& paramcurve ,
		int cellDim, const QuadratureFamily& quad) : StdFwkEvalMediator(mesh, cellDim) ,
    numEvaluationCases_(-1),
    quad_(quad),
    numQuadPtsForMaxCell_(0) ,
    paramcurve_(paramcurve) {

	// set the cell types
	maxCellType_ = mesh.cellType(mesh.spatialDim() );
	curveCellType_ = mesh.cellType(mesh.spatialDim()-1);

	// detect the (constant) quadrature points (we must have the same nr quad points per cell!!!)
	Array<Point> quadPoints;
	Array<double> quadWeights;
	quad.getPoints(curveCellType_ , quadPoints , quadWeights );
	numQuadPtsForMaxCell_ = quadWeights.size();
}
Ejemplo n.º 6
0
int main(int argc, char** argv)
{
  
  try
  {
    GlobalMPISession session(&argc, &argv);

    TimeMonitor t(totalTimer());

    int pMax = 1;
    int dim=2;

    CellType cellType = TriangleCell;

    Point a = Point(0.0, 0.0);
    Point b = Point(1.0, 0.0);
    Point c = Point(0.0, 1.0);
    CellJacobianBatch JBatch;
    JBatch.resize(1, 2, 2);
    double* J = JBatch.jVals(0);
    J[0] = b[0] - a[0];
    J[1] = c[0] - a[0];
    J[2] = b[1] - a[1];
    J[3] = c[1] - a[1];


    bool isInternalBdry=false;

    /* ------ evaluate Lagrange and FIAT-Lagrange at the vertices */
    Array<Point> verts = tuple(a,b,c);
    BasisFamily lagrange = new Lagrange(1);
    BasisFamily fiatLagrange = new Lagrange(1);
      
    MultiIndex d0(0,0,0);
    MultiIndex dx(1,0,0);
    MultiIndex dy(0,1,0);

    Array<Array<Array<double> > > result;

    Array<int> dummy;

    std::cerr << "------ Evaluating bases at vertices ----------" << std::endl
         << std::endl;

    std::cerr << "Evaluating phi(vert) with FIAT-Lagrange" << std::endl;
    fiatLagrange.ptr()->refEval(cellType, verts, d0, result);
    std::cerr << "results = " << result << std::endl << std::endl;

    std::cerr << "Evaluating phi(vert) with Lagrange" << std::endl;
    lagrange.ptr()->refEval(cellType, verts, d0, result);
    std::cerr << "results = " << result << std::endl << std::endl;

    std::cerr << std::endl ;

    std::cerr << "Evaluating Dx*phi(vert) with FIAT-Lagrange" << std::endl;
    fiatLagrange.ptr()->refEval(cellType, verts, dx, result);
    std::cerr << "results = " << result << std::endl << std::endl;

    std::cerr << "Evaluating Dx*phi(vert) with Lagrange" << std::endl;
    lagrange.ptr()->refEval(cellType, verts, dx, result);
    std::cerr << "results = " << result << std::endl << std::endl;

    std::cerr << std::endl ;
      
    std::cerr << "Evaluating Dy*phi(vert) with FIAT-Lagrange" << std::endl;
    fiatLagrange.ptr()->refEval(cellType, verts, dy, result);
    std::cerr << "results = " << result << std::endl << std::endl;

    std::cerr << "Evaluating Dy*phi(vert) with Lagrange" << std::endl;
    lagrange.ptr()->refEval(cellType, verts, dy, result);
    std::cerr << "results = " << result << std::endl << std::endl;

      

    /* --------- evaluate integrals over elements ----------- */
      
    RCP<Array<double> > A = rcp(new Array<double>());
          
    QuadratureFamily quad = new GaussianQuadrature(4);
    Array<double> quadWeights;
    Array<Point> quadPts;
    quad.getPoints(cellType, quadPts, quadWeights);
    int nQuad = quadPts.size();

    Array<double> coeff(nQuad);
    for (int i=0; i<nQuad; i++) 
    {
      double s = quadPts[i][0];
      double t = quadPts[i][1];
      double x = a[0] + J[0]*s + J[1]*t;
      double y = a[1] + J[2]*s + J[3]*t;
      coeff[i] = x*y;
    }
    const double* const f = &(coeff[0]);

    std::cerr << std::endl << std::endl 
         << "---------------- One-forms --------------------" 
         << std::endl << std::endl;
    for (int p=1; p<=pMax; p++)
    {
      BasisFamily P = new Lagrange(p);
      for (int dp=0; dp<=1; dp++)
      {
        if (dp > p) continue;
        Tabs tab0;
        std::cerr << tab0 << "test function deriv order = " << dp << std::endl;
        int numTestDir = 1;
        if (dp==1) numTestDir = dim;
        for (int t=0; t<numTestDir; t++)
        {
          int alpha = t;
          Tabs tab;
          QuadratureIntegral ref(dim, cellType, dim, cellType, P, alpha, dp, quad, isInternalBdry);
          A->resize(ref.nNodesTest());
          ref.transformOneForm(JBatch, JBatch, dummy, f, A);
          std::cerr << tab << "test deriv direction =" << t << std::endl;
          std::cerr << tab << "transformed local vector: " << std::endl;
          std::cerr << tab << "{";
          for (int r=0; r<ref.nNodesTest(); r++)
          {
            if (r!=0) std::cerr << ", ";
            std::cerr << (*A)[r];
          }
          std::cerr << "}" << std::endl << std::endl;
        }
      }
    }

    std::cerr << std::endl << std::endl 
         << "---------------- Two-forms --------------------" 
         << std::endl << std::endl;
    for (int p=1; p<=pMax; p++)
    {
      BasisFamily P = new Lagrange(p);
      for (int q=1; q<=pMax; q++)
      {
        BasisFamily Q = new Lagrange(q);
        for (int dp=0; dp<=1; dp++)
        {
          if (dp > p) continue;
          Tabs tab0;
          std::cerr << tab0 << "test function deriv order = " << dp << std::endl;
          for (int dq=0; dq<=1; dq++)
          {
            if (dq > q) continue;
            Tabs tab1;
            std::cerr << tab1 
                 << "unk function deriv order = " << dq << std::endl;
            int numTestDir = 1;
            if (dp==1) numTestDir = dim;
            for (int t=0; t<numTestDir; t++)
            {
              int alpha = t;
              int numUnkDir = 1;
              if (dq==1) numUnkDir = dim;
              for (int u=0; u<numUnkDir; u++)
              {
                Tabs tab;
                int beta = u;
                QuadratureIntegral ref(dim, cellType, dim, cellType, P, alpha, 
                  dp, Q, beta, dq, quadd, isInternalBdry);
                A->resize(ref.nNodesTest()*ref.nNodesUnk());
                ref.transformTwoForm(JBatch, JBatch, dummy, f, A);

                std::cerr << tab << "test deriv direction =" << 
                  t << ", unk deriv direction =" << u << std::endl;
                std::cerr << tab << "transformed local stiffness matrix" << std::endl;
                std::cerr << tab << "{";

                for (int r=0; r<ref.nNodesTest(); r++)
                {
                  if (r!=0) std::cerr << ", ";
                  std::cerr << "{";
                  for (int c=0; c<ref.nNodesUnk(); c++)
                  {
                    if (c!=0) std::cerr << ", ";
                    std::cerr << chop((*A)[r + ref.nNodesTest()*c]);
                  }
                  std::cerr << "}";
                }
                std::cerr << "}" << std::endl << std::endl;
              }
            }
          }
        }
      }
    }
    TimeMonitor::summarize();

  }
	catch(std::exception& e)
  {
    std::cerr << e.what() << std::endl;
  }
}
void checkbasis( BasisFamily &b1 , BasisFamily &b2 )
{
  int maxDim=3;
  double tol = 1.0e-13;
  int maxDiffOrder = 0;
  int numErrors = 0;
  QuadratureFamily quad = new GaussianQuadrature(4);
  
  for (int spatialDim=1; spatialDim<=maxDim; spatialDim++) {
    std::cerr << "\t" << "spatial dimension =" << spatialDim << std::endl;
    for (int cellDim=0; cellDim<=spatialDim; cellDim++) { 
      std::cerr << "\t\t" << "cell dimension =" << cellDim << std::endl;
      CellType cellType;
      if (cellDim==0) cellType=PointCell;
      if (cellDim==1) cellType=LineCell;
      if (cellDim==2) cellType=TriangleCell;
      if (cellDim==3) cellType=TetCell;
      
      Array<Point> qPts;
      Array<double> qWts;
      quad.getPoints(cellType, qPts, qWts);
      
      for (int d=0; d<=maxDiffOrder; d++) {
	if (cellDim==0 && d>0) continue;
	cerr << "\t\t\t" << "differentiation order = " << d << std::endl;
	for (int dir=0; dir<iPow(cellDim, d); dir++) {
	  std::cerr << "\t\t\t\t" << "direction = " << dir << std::endl;
	  MultiIndex mi;
	  mi[dir]=d;
	  Array<Array<double> > values1;
	  Array<Array<double> > values2;
	  std::cerr << "\t\t\t\t" << "computing basis1...";
	  b1.ptr()->refEval(spatialDim, cellType, qPts, mi, values1);
	  std::cerr << "done" << std::endl << "\t\t\t\t" << "computing basis2...";
	  b2.ptr()->refEval(spatialDim, cellType, qPts, mi, values2);
	  std::cerr << "done" << std::endl;
	  int nNodes1 = b1.ptr()->nNodes(spatialDim, cellType);
	  int nNodes2 = b2.ptr()->nNodes(spatialDim, cellType);
	  std::cerr << "\t\t\t\t" << "num nodes: basis1=" << nNodes1
	       << " basis2=" << nNodes2 << std::endl;
	  if (nNodes1 != nNodes2) {	
	    std::cerr << "******** ERROR: node counts should be equal" << std::endl;
	    numErrors++;
	    continue;
	  }
	  if (values1.size() != values2.size()) {
	    std::cerr << "******** ERROR: value array outer sizes should be equal" << std::endl;
	    numErrors++;
	    continue;
	  }
	  if (values1.size() != qPts.size()) {
	    std::cerr << "******** ERROR: value array outer size should be equal to number of quad points" << std::endl;
	    numErrors++;
	    continue;
	  }
	  for (int q=0; q<qPts.length(); q++) {
	    if (values1[q].length() != nNodes1) {
	      std::cerr << "******** ERROR: value array inner size should be equal to number of nodes" << std::endl;
	      numErrors++;
	      continue;
	    }
	    std::cerr << "\t\t\t\t\t" << "quad point q=" << q << " pt=" << qPts[q]
		 << std::endl;
	    for (int n=0; n<nNodes1; n++) {
	      std::cerr << "\t\t\t\t\t\t" << "node n=" << n << " phi1="
		   << values1[q][n]
		   << " phi2=" << values2[q][n]
		   << " |phi1-phi2|=" << fabs(values1[q][n]-values2[q][n])
		   << std::endl;
	      if (fabs(values1[q][n]-values2[q][n]) > tol) {
		cout << "ERROR" << std::endl; numErrors++;
	      }
	    }
	  }
	}
      }
    }
  }    
  std::cerr << std::endl << std::endl << "Summary: detected " << numErrors << " errors " << std::endl;
}
MaximalQuadratureIntegral::MaximalQuadratureIntegral(
  const CellType& cellType,
  const BasisFamily& testBasis,
  int alpha,
  int testDerivOrder,
  const BasisFamily& unkBasis,
  int beta,
  int unkDerivOrder,
  const QuadratureFamily& quad,
  const ParametrizedCurve& globalCurve,
  const Mesh& mesh,
  int verb)
  : ElementIntegral(dimension(cellType), cellType, dimension(cellType),
    cellType, testBasis,
    alpha, testDerivOrder, unkBasis, beta, unkDerivOrder, true,
    globalCurve, mesh, 
    verb),
    quad_(quad),
    quadPts_(),
    quadWeights_(),
    W_(),
    useSumFirstMethod_(true)
{
  Tabs tab0(0);
  
  SUNDANCE_MSG1(setupVerb(), tab0 << "MaximalQuadratureIntegral ctor for 2-form");
  if (setupVerb()) describe(Out::os());
  assertBilinearForm();

  // store the quadrature points and weights  
  quad.getPoints(cellType, quadPts_, quadWeights_);
  int nQuad = quadPts_.size();

  W_.resize(nQuad * nRefDerivTest() * nNodesTest()  
    * nRefDerivUnk() * nNodesUnk());


  /* compute the basis functions */
  Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());
  Array<Array<Array<Array<double> > > > unkBasisVals(nRefDerivUnk());
  

  for (int r=0; r<nRefDerivTest(); r++)
  {
    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++)
  {
    unkBasisVals[r].resize(unkBasis.dim());
    MultiIndex mi;
    if (unkDerivOrder==1) mi[r] = 1;
    SpatialDerivSpecifier deriv(mi);
    unkBasis.refEval(evalCellType(), 
      quadPts_, deriv, unkBasisVals[r], setupVerb());
  }
  

  int vecComp = 0;
  /* form the products of basis functions at each quad pt */
  W_ACI_F2_.resize(nQuad);
  for (int q=0; q<nQuad; q++)
  {
    W_ACI_F2_[q].resize(nRefDerivTest());
    for (int t=0; t<nRefDerivTest(); t++)
    {
      W_ACI_F2_[q][t].resize(nNodesTest());
      for (int nt=0; nt<nNodesTest(); nt++)
      {
        W_ACI_F2_[q][t][nt].resize(nRefDerivUnk());
        for (int u=0; u<nRefDerivUnk(); u++)
        {
          W_ACI_F2_[q][t][nt][u].resize(nNodesUnk());
          for (int nu=0; nu<nNodesUnk(); nu++)
          {
            wValue(q, t, nt, u, nu)
              = chop(quadWeights_[q] * testBasisVals[t][vecComp][q][nt] 
                * unkBasisVals[u][vecComp][q][nu]);
            W_ACI_F2_[q][t][nt][u][nu] =
            	chop(testBasisVals[t][vecComp][q][nt] * unkBasisVals[u][vecComp][q][nu]);
          }
        }
      }
    }
  }

  addFlops(3*nQuad*nRefDerivTest()*nNodesTest()*nRefDerivUnk()*nNodesUnk()
    + W_.size());
  for (int i=0; i<W_.size(); i++) W_[i] = chop(W_[i]);

}
MaximalQuadratureIntegral::MaximalQuadratureIntegral(
  const CellType& cellType,
  const BasisFamily& testBasis,
  int alpha,
  int testDerivOrder,
  const QuadratureFamily& quad,
  const ParametrizedCurve& globalCurve,
  const Mesh& mesh,
  int verb)
  : ElementIntegral(dimension(cellType), cellType, dimension(cellType),
    cellType, testBasis,
    alpha, testDerivOrder, true, globalCurve, mesh, verb),
    quad_(quad),
    quadPts_(),
    quadWeights_(),
    W_(),
    useSumFirstMethod_(true)
{
  Tabs tab0(0);
  
  SUNDANCE_MSG1(setupVerb(), tab0 << "MaximalQuadratureIntegral ctor for 1-form");
  if (setupVerb()) describe(Out::os());
  assertLinearForm();

  
  SUNDANCE_MSG1(setupVerb(), tab0 << "quadrature family=" << quad);  
  
  quad.getPoints(cellType, quadPts_, quadWeights_);
  int nQuad = quadPts_.size();
      
  W_.resize(nQuad * nRefDerivTest() * nNodesTest());

  SUNDANCE_MSG1(setupVerb(), tab0 << "num nodes for test function " << nNodesTest());

  Array<Array<Array<Array<double> > > > testBasisVals(nRefDerivTest());

  for (int r=0; r<nRefDerivTest(); r++)
  {
    Tabs tab3;
    SUNDANCE_MSG1(setupVerb(), tab3 
      << "evaluating basis functions for ref deriv direction " << r);
    MultiIndex mi;
    testBasisVals[r].resize(testBasis.dim());
    if (testDerivOrder==1) mi[r] = 1;
    SpatialDerivSpecifier deriv(mi);
    testBasis.refEval(evalCellType(), quadPts_, deriv, 
      testBasisVals[r], setupVerb());
  }

  int vecComp = 0;
  W_ACI_F1_.resize(nQuad);
  for (int q=0; q<nQuad; q++)
  {
    W_ACI_F1_[q].resize(nRefDerivTest());
    for (int t=0; t<nRefDerivTest(); t++)
    {
      W_ACI_F1_[q][t].resize(nNodesTest());
      for (int nt=0; nt<nNodesTest(); nt++)
      {
        wValue(q, t, nt) 
          = chop(quadWeights_[q] * testBasisVals[t][vecComp][q][nt]) ;
        W_ACI_F1_[q][t][nt] = chop(testBasisVals[t][vecComp][q][nt]);
      }
    }
  }

  addFlops(2*nQuad*nRefDerivTest()*nNodesTest());
}
Ejemplo n.º 10
0
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");
}
Ejemplo n.º 11
0
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");
}