void CurveEvalMediator
::evalDiscreteFuncElement(const DiscreteFuncElement* expr,
  const Array<MultiIndex>& multiIndices,
  Array<RCP<EvalVector> >& vec) const
{

	  int verbo = dfVerb();
	  Tabs tab1;
	  SUNDANCE_MSG2(verbo , tab1
	    << "CurveEvalMediator evaluating Discrete Function expr " << expr->toString());

	  const DiscreteFunctionData* f = DiscreteFunctionData::getData(expr);

	  TEUCHOS_TEST_FOR_EXCEPTION(f==0, std::logic_error,
	    "QuadratureEvalMediator::evalDiscreteFuncElement() called "
	    "with expr that is not a discrete function");

	  SUNDANCE_MSG2(verbo , tab1 << "After casting DiscreteFunctionData" << expr->toString());

	  RCP<Array<Array<double> > > localValues;
	  Array<int> cellLIDs_tmp(1);
	  Array<Point> phyPoints;
	  Array<Point> refPoints;
	  Array<Point> refDevs;
	  Array<Point> refNormal;
	  int nCells = cellLID()->size();

      RCP<const MapStructure> mapStruct;
	  int myIndex = expr->myIndex();
	  int nQuad = numQuadPtsForMaxCell_;
	  Array<int> k(multiIndices.size(),0);

	  Teuchos::BLAS<int,double> blas;

	  SUNDANCE_MSG2(verbo , tab1 << "After declaring BLAS: " << expr->toString());

	  // resize correctly the result vector
	  for (int i=0; i<multiIndices.size(); i++)
	  {
		  vec[i]->resize(nCells*nQuad);
	  }

	  // loop over each cell
	  for (int c=0; c<nCells; c++)
	  {
			int maxCellLID = (*cellLID())[c];
		    localValues = rcp(new Array<Array<double> >());
			SUNDANCE_MSG2(verbo , tab1 <<  "Cell:" << c <<  " of " << nCells << " , maxCellLID:" << maxCellLID );

			cellLIDs_tmp[0] = (*cellLID())[c];
			SUNDANCE_MSG2(verbo , tab1 <<  " Before calling f->getLocalValues:" << cellLIDs_tmp.size() <<
					" f==0 : " << (f==0) << " tmp:" << f->mesh().spatialDim());
			// - get local values from the DiscreteFunctionElementData
			mapStruct = f->getLocalValues(maxCellDim(), cellLIDs_tmp , *localValues);
			SUNDANCE_MSG2(verbo , tab1 <<  " After getting mapStruct:" << maxCellLID );
			SUNDANCE_MSG2(verbo , tab1 <<  " mapStruct->numBasisChunks():" << mapStruct->numBasisChunks() );
		    int chunk = mapStruct->chunkForFuncID(myIndex);
		    int funcIndex = mapStruct->indexForFuncID(myIndex);
		    int nFuncs = mapStruct->numFuncs(chunk);
			SUNDANCE_MSG2(verbo , tab1 <<  " chunk:" << chunk );
			SUNDANCE_MSG2(verbo , tab1 <<  " funcIndex:" << funcIndex );
			SUNDANCE_MSG2(verbo , tab1 <<  " nFuncs:" << nFuncs );

			// the chunk of the function
		    BasisFamily basis = rcp_dynamic_cast<BasisFamilyBase>(mapStruct->basis(chunk));
		    int nNodesTotal = basis.nReferenceDOFsWithFacets(maxCellType(), maxCellType());

			// - get intersection (reference)points from the mesh (if not existent than compute them)
			if ( mesh().hasCurvePoints( maxCellLID , paramcurve_.myID() ))
			{
				mesh().getCurvePoints( maxCellLID , paramcurve_.myID() , refPoints , refDevs , refNormal );
			}
			else // we have to calculate now the points
			{
				// calculate the intersection points
				CurveIntegralCalc::getCurveQuadPoints( maxCellType_ , maxCellLID , mesh() , paramcurve_ , quad_ ,
						refPoints, refDevs , refNormal);
				// store the intersection point in the mesh
				mesh().setCurvePoints( maxCellLID , paramcurve_.myID() , refPoints, refDevs , refNormal );
			}

			// loop over each multi-index
			SUNDANCE_MSG2(verbo , tab1 << " multiIndices.size()" << multiIndices.size() );
			for (int i=0; i<multiIndices.size(); i++)
			{

				int nDerivResults = 1;
				if ( multiIndices[i].order() == 1 ) nDerivResults = maxCellDim();

			    int pDir = 0;
			    int derivNum = 1;
				MultiIndex mi;
				SUNDANCE_MSG2(verbo , tab1 << " before asking anything i = " << i);
				SUNDANCE_MSG2(verbo , tab1 << " multiindex order : " << multiIndices[i].order());
				SUNDANCE_MSG2(verbo , tab1 << " multiindex : " << multiIndices[i] );
				if (multiIndices[i].order() > 0){
					pDir = multiIndices[i].firstOrderDirection();
					mi[pDir] = 1;
					derivNum = mesh().spatialDim();
				}

				Array<Array<double> > result(nQuad*derivNum);
				Array<Array<Array<double> > > tmp;

				int offs = nNodesTotal;

				// resize the result vector
				for (int deriv = 0 ; deriv < derivNum ; deriv++)
				{
                    // test weather we have to compute derivative
					if (multiIndices[i].order() > 0){
						// in case of derivatives we set one dimension
						MultiIndex mi_tmp;
						mi_tmp[deriv] = 1;
						SpatialDerivSpecifier deriv(mi_tmp);
						SUNDANCE_MSG2(verbo , tab1 << "computing derivatives : " << deriv << " on reference cell ");
						basis.refEval( maxCellType_ , refPoints , deriv, tmp , verbo );
					} else {
						SpatialDerivSpecifier deriv(mi);
						 // --- end eval basis functions
						SUNDANCE_MSG2(verbo , tab1 << "computing values reference cell ");
						basis.refEval( maxCellType_ , refPoints , deriv, tmp , verbo );
					}

					SUNDANCE_MSG2(verbo , tab1 << "resize result vector , offs:" << offs);
					for (int q=0; q<nQuad; q++){
						result[nQuad*deriv + q].resize(offs);
					}

				   	// copy the result in an other format
					SUNDANCE_MSG2(verbo , tab1 << "copy results ");
					int  offs1 = 0;
				    for (int q=0; q<nQuad; q++)
				    {
				    	offs1 = 0;
						for (int d=0; d<basis.dim(); d++)
						{
						   int nNodes = tmp[d][q].size();
					       for (int n=0; n<nNodes; n++ , offs1++ )
					       {
					          	result[nQuad*deriv + q][offs1] = tmp[d][q][n];
					       }
					    }
					}
				}// loop over all dimensional derivative

                // multiply the local results with the coefficients, (matrix vector OP)
				SUNDANCE_MSG2(verbo , tab1 << "summing up values , funcIndex:" << funcIndex << " offs:" << offs);
				for (int deriv = 0 ; deriv < derivNum ; deriv++)
				{
					for (int q=0; q<nQuad; q++)
					{
						double sum = 0.0;
						// sum over nodes
						for (int n = 0 ; n < offs ; n++){
							sum = sum + result[nQuad*deriv + q][n] * (*localValues)[chunk][funcIndex*offs + n];
						}
						// sum up the result in the 0th element
						result[nQuad*deriv + q][0] = sum;
					}
				}

			    // multiply the result if necesary with the inverse of the Jacobian
				const CellJacobianBatch& J = JTrans();
				if (mi.order()==1)
				{
					Tabs tab1;
					Tabs tab2;
				    SUNDANCE_MSG2(verbo, tab2 << "Jacobian batch nCells=" << J.numCells());
			        SUNDANCE_MSG2(verbo, tab2 << "Jacobian batch cell dim=" << J.cellDim());
			        SUNDANCE_MSG2(verbo, tab2 << "Jacobian batch spatial dim=" << J.spatialDim());
                    // we just multiply the derivative direction component
			        Array<double> invJ;
			        J.getInvJ(c, invJ);
					for (int q=0; q<nQuad; q++)
					{
						double sum = 0.0;
						for (int deriv = 0 ; deriv < derivNum ; deriv++)
						{
							// multiply one row from the J^{-T} matrix with the gradient vector
							sum = sum + result[nQuad*deriv + q][0] * invJ[derivNum*pDir + deriv];
						}
						// the resulting derivative on the physical cell in the "pDir" direction
						result[q][0] = sum;
					}
			    }

				// --- just copy the result to the "vec" back, the result should be in the  "result[q][0]" place----
				//SUNDANCE_MSG2(verbo , tab1 << "copy results back ");
				double* vecPtr = vec[i]->start();
				for (int q=0; q<nQuad; q++, k[i]++)
				{
					vecPtr[k[i]] = result[q][0];
				}
				SUNDANCE_MSG2(verbo , tab1 << " END copy results back ");
			} // --- end loop multiindex
			SUNDANCE_MSG2(verbo , tab1 << " END loop over multiindex ");
	}// --- end loop over cells
	SUNDANCE_MSG2(verbo , tab1 << " END loop over cells ");
}
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]);

}
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");
}
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());
}
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");
}