void values(FieldContainer<double> &values, BasisCachePtr sliceBasisCache) { vector<GlobalIndexType> sliceCellIDs = sliceBasisCache->cellIDs(); Teuchos::Array<int> dim; values.dimensions(dim); dim[0] = 1; // one cell Teuchos::Array<int> offset(dim.size()); for (int cellOrdinal = 0; cellOrdinal < sliceCellIDs.size(); cellOrdinal++) { offset[0] = cellOrdinal; int enumeration = values.getEnumeration(offset); FieldContainer<double>valuesForCell(dim,&values[enumeration]); GlobalIndexType sliceCellID = sliceCellIDs[cellOrdinal]; int numPoints = sliceBasisCache->getPhysicalCubaturePoints().dimension(1); int spaceDim = sliceBasisCache->getPhysicalCubaturePoints().dimension(2); FieldContainer<double> spaceTimePhysicalPoints(1,numPoints,spaceDim+1); for (int ptOrdinal=0; ptOrdinal<numPoints; ptOrdinal++) { for (int d=0; d<spaceDim; d++) { spaceTimePhysicalPoints(0,ptOrdinal,d) = sliceBasisCache->getPhysicalCubaturePoints()(cellOrdinal,ptOrdinal,d); } spaceTimePhysicalPoints(0,ptOrdinal,spaceDim) = _t; } GlobalIndexType cellID = _cellIDMap[sliceCellID]; BasisCachePtr spaceTimeBasisCache = BasisCache::basisCacheForCell(_spaceTimeMesh, cellID); FieldContainer<double> spaceTimeRefPoints(1,numPoints,spaceDim+1); CamelliaCellTools::mapToReferenceFrame(spaceTimeRefPoints, spaceTimePhysicalPoints, _spaceTimeMesh, cellID); spaceTimeRefPoints.resize(numPoints,spaceDim+1); spaceTimeBasisCache->setRefCellPoints(spaceTimeRefPoints); _spaceTimeFunction->values(valuesForCell, spaceTimeBasisCache); } }
void BilinearFormUtility::weightCellBasisValues(FieldContainer<double> &basisValues, const FieldContainer<double> &weights, int offset) { // weights are (numCells, offset+numFields) // basisValues are (numCells, numFields, ...) int numCells = basisValues.dimension(0); int numFields = basisValues.dimension(1); Teuchos::Array<int> dimensions; basisValues.dimensions(dimensions); int numAffectedValues = 1; for (int dimIndex=2; dimIndex<dimensions.size(); dimIndex++) { numAffectedValues *= dimensions[dimIndex]; } Teuchos::Array<int> index(dimensions.size(),0); for (int cellIndex=0; cellIndex < numCells; cellIndex++) { index[0] = cellIndex; for (int fieldIndex=0; fieldIndex < numFields; fieldIndex++) { index[1] = fieldIndex; int enumIndex = basisValues.getEnumeration(index); for (int valIndex=enumIndex; valIndex < numAffectedValues + enumIndex; valIndex++) { basisValues[valIndex] *= weights(cellIndex,fieldIndex+offset); } } } }
void SpatiallyFilteredFunction<Scalar>::values(FieldContainer<Scalar> &values, BasisCachePtr basisCache) { // cout << "Entered SpatiallyFilteredFunction<Scalar>::values()\n"; int numCells = values.dimension(0); int numPoints = values.dimension(1); values.initialize(0.0); Teuchos::Array<int> dim; values.dimensions(dim); Teuchos::Array<int> fValuesDim = dim; int entriesPerPoint = 1; for (int d=2; d<values.rank(); d++) { entriesPerPoint *= dim[d]; dim[d] = 0; // clear so that these indices point to the start of storage for (cellIndex,ptIndex) } FieldContainer<bool> pointsMatch(numCells,numPoints); if (_sf->matchesPoints(pointsMatch,basisCache)) // SOME point matches { // cout << "pointsMatch:\n" << pointsMatch; FieldContainer<Scalar> fValues(fValuesDim); _f->values(fValues,basisCache); for (int cellIndex=0; cellIndex<numCells; cellIndex++) { dim[0] = cellIndex; for (int ptIndex=0; ptIndex<numPoints; ptIndex++) { dim[1] = ptIndex; if (pointsMatch(cellIndex,ptIndex)) { Scalar* value = &values[values.getEnumeration(dim)]; Scalar* fValue = &fValues[fValues.getEnumeration(dim)]; for (int entryIndex=0; entryIndex<entriesPerPoint; entryIndex++) { *value++ = *fValue++; } } } } } }
void MPIWrapper::entryWiseSum(FieldContainer<GlobalIndexType> &values) { #ifdef HAVE_MPI // cast to long long: Teuchos::Array<int> dim; values.dimensions(dim); FieldContainer<long long> valuesLongLong(dim); for (int i=0; i<values.size(); i++) { valuesLongLong[i] = (long long) values[i]; } Epetra_MpiComm Comm(MPI_COMM_WORLD); FieldContainer<long long> valuesLongLongCopy = valuesLongLong; // it appears this copy is necessary Comm.SumAll(&valuesLongLongCopy[0], &valuesLongLong[0], valuesLongLong.size()); // copy back to original container: for (int i=0; i<values.size(); i++) { values[i] = (GlobalIndexType) valuesLongLong[i]; } #else #endif }
int main(int argc, char *argv[]) { #ifdef HAVE_MPI Teuchos::GlobalMPISession mpiSession(&argc, &argv,0); int rank=mpiSession.getRank(); int numProcs=mpiSession.getNProc(); #else int rank = 0; int numProcs = 1; #endif int polyOrder = 3; int pToAdd = 2; // for tests // define our manufactured solution or problem bilinear form: double epsilon = 1e-2; bool useTriangles = false; FieldContainer<double> quadPoints(4,2); quadPoints(0,0) = 0.0; // x1 quadPoints(0,1) = 0.0; // y1 quadPoints(1,0) = 1.0; quadPoints(1,1) = 0.0; quadPoints(2,0) = 1.0; quadPoints(2,1) = 1.0; quadPoints(3,0) = 0.0; quadPoints(3,1) = 1.0; int H1Order = polyOrder + 1; int horizontalCells = 1, verticalCells = 1; double energyThreshold = 0.2; // for mesh refinements double nonlinearStepSize = 0.5; double nonlinearRelativeEnergyTolerance = 0.015; // used to determine convergence of the nonlinear solution //////////////////////////////////////////////////////////////////// // SET UP PROBLEM //////////////////////////////////////////////////////////////////// Teuchos::RCP<BurgersBilinearForm> oldBurgersBF = Teuchos::rcp(new BurgersBilinearForm(epsilon)); // new-style bilinear form definition VarFactory varFactory; VarPtr uhat = varFactory.traceVar("\\widehat{u}"); VarPtr beta_n_u_minus_sigma_hat = varFactory.fluxVar("\\widehat{\\beta_n u - \\sigma_n}"); VarPtr u = varFactory.fieldVar("u"); VarPtr sigma1 = varFactory.fieldVar("\\sigma_1"); VarPtr sigma2 = varFactory.fieldVar("\\sigma_2"); VarPtr tau = varFactory.testVar("\\tau",HDIV); VarPtr v = varFactory.testVar("v",HGRAD); BFPtr bf = Teuchos::rcp( new BF(varFactory) ); // create a pointer to a new mesh: Teuchos::RCP<Mesh> mesh = Mesh::buildQuadMesh(quadPoints, horizontalCells, verticalCells, bf, H1Order, H1Order+pToAdd, useTriangles); mesh->setPartitionPolicy(Teuchos::rcp(new ZoltanMeshPartitionPolicy("HSFC"))); Teuchos::RCP<Solution> backgroundFlow = Teuchos::rcp(new Solution(mesh, Teuchos::rcp((BC*)NULL) , Teuchos::rcp((RHS*)NULL), Teuchos::rcp((DPGInnerProduct*)NULL))); // create null solution oldBurgersBF->setBackgroundFlow(backgroundFlow); // tau parts: // 1/eps (sigma, tau)_K + (u, div tau)_K - (u_hat, tau_n)_dK bf->addTerm(sigma1 / epsilon, tau->x()); bf->addTerm(sigma2 / epsilon, tau->y()); bf->addTerm(u, tau->div()); bf->addTerm( - uhat, tau->dot_normal() ); vector<double> e1(2); // (1,0) e1[0] = 1; vector<double> e2(2); // (0,1) e2[1] = 1; FunctionPtr u_prev = Teuchos::rcp( new PreviousSolutionFunction(backgroundFlow, u) ); FunctionPtr beta = e1 * u_prev + Teuchos::rcp( new ConstantVectorFunction( e2 ) ); // v: // (sigma, grad v)_K - (sigma_hat_n, v)_dK - (u, beta dot grad v) + (u_hat * n dot beta, v)_dK bf->addTerm( sigma1, v->dx() ); bf->addTerm( sigma2, v->dy() ); bf->addTerm( -u, beta * v->grad()); bf->addTerm( beta_n_u_minus_sigma_hat, v); // ==================== SET INITIAL GUESS ========================== mesh->registerSolution(backgroundFlow); map<int, Teuchos::RCP<AbstractFunction> > functionMap; functionMap[BurgersBilinearForm::U] = Teuchos::rcp(new InitialGuess()); functionMap[BurgersBilinearForm::SIGMA_1] = Teuchos::rcp(new ZeroFunction()); functionMap[BurgersBilinearForm::SIGMA_2] = Teuchos::rcp(new ZeroFunction()); backgroundFlow->projectOntoMesh(functionMap); // ==================== END SET INITIAL GUESS ========================== // compare stiffness matrices for first linear step: int trialOrder = 1; pToAdd = 0; int testOrder = trialOrder + pToAdd; CellTopoPtr quadTopoPtr = Teuchos::rcp(new shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() )); DofOrderingFactory dofOrderingFactory(bf); DofOrderingPtr testOrdering = dofOrderingFactory.testOrdering(testOrder, *quadTopoPtr); DofOrderingPtr trialOrdering = dofOrderingFactory.trialOrdering(trialOrder, *quadTopoPtr); int numCells = 1; // just use testOrdering for both trial and test spaces (we only use to define BasisCache) ElementTypePtr elemType = Teuchos::rcp( new ElementType(trialOrdering, testOrdering, quadTopoPtr) ); BasisCachePtr basisCache = Teuchos::rcp( new BasisCache(elemType) ); quadPoints.resize(1,quadPoints.dimension(0),quadPoints.dimension(1)); basisCache->setPhysicalCellNodes(quadPoints,vector<int>(1),true); // true: do create side cache FieldContainer<double> cellSideParities(numCells,quadTopoPtr->getSideCount()); cellSideParities.initialize(1.0); // not worried here about neighbors actually having opposite parity -- just want the two BF implementations to agree... FieldContainer<double> expectedValues(numCells, testOrdering->totalDofs(), trialOrdering->totalDofs() ); FieldContainer<double> actualValues(numCells, testOrdering->totalDofs(), trialOrdering->totalDofs() ); oldBurgersBF->stiffnessMatrix(expectedValues, elemType, cellSideParities, basisCache); bf->stiffnessMatrix(actualValues, elemType, cellSideParities, basisCache); // compare beta's as well FieldContainer<double> expectedBeta = oldBurgersBF->getBeta(basisCache); Teuchos::Array<int> dim; expectedBeta.dimensions(dim); FieldContainer<double> actualBeta(dim); beta->values(actualBeta,basisCache); double tol = 1e-14; double maxDiff; if (rank == 0) { if ( ! TestSuite::fcsAgree(expectedBeta,actualBeta,tol,maxDiff) ) { cout << "Test failed: old Burgers beta differs from new; maxDiff " << maxDiff << ".\n"; cout << "Old beta: \n" << expectedBeta; cout << "New beta: \n" << actualBeta; } else { cout << "Old and new Burgers beta agree!!\n"; } if ( ! TestSuite::fcsAgree(expectedValues,actualValues,tol,maxDiff) ) { cout << "Test failed: old Burgers stiffness differs from new; maxDiff " << maxDiff << ".\n"; cout << "Old: \n" << expectedValues; cout << "New: \n" << actualValues; cout << "TrialDofOrdering: \n" << *trialOrdering; cout << "TestDofOrdering:\n" << *testOrdering; } else { cout << "Old and new Burgers stiffness agree!!\n"; } } // define our inner product: // Teuchos::RCP<BurgersInnerProduct> ip = Teuchos::rcp( new BurgersInnerProduct( bf, mesh ) ); // function to scale the squared guy by epsilon/h FunctionPtr epsilonOverHScaling = Teuchos::rcp( new EpsilonScaling(epsilon) ); IPPtr ip = Teuchos::rcp( new IP ); ip->addTerm(tau); ip->addTerm(tau->div()); ip->addTerm( epsilonOverHScaling * v ); ip->addTerm( sqrt(sqrt(epsilon)) * v->grad() ); ip->addTerm( beta * v->grad() ); // use old IP instead, for now... Teuchos::RCP<BurgersInnerProduct> oldIP = Teuchos::rcp( new BurgersInnerProduct( oldBurgersBF, mesh ) ); expectedValues.resize(numCells, testOrdering->totalDofs(), testOrdering->totalDofs() ); actualValues.resize (numCells, testOrdering->totalDofs(), testOrdering->totalDofs() ); BasisCachePtr ipBasisCache = Teuchos::rcp( new BasisCache(elemType, true) ); // true: test vs. test ipBasisCache->setPhysicalCellNodes(quadPoints,vector<int>(1),false); // false: don't create side cache oldIP->computeInnerProductMatrix(expectedValues,testOrdering,ipBasisCache); ip->computeInnerProductMatrix(actualValues,testOrdering,ipBasisCache); tol = 1e-14; maxDiff = 0.0; if (rank==0) { if ( ! TestSuite::fcsAgree(expectedValues,actualValues,tol,maxDiff) ) { cout << "Test failed: old inner product differs from new IP; maxDiff " << maxDiff << ".\n"; cout << "Old: \n" << expectedValues; cout << "New IP: \n" << actualValues; cout << "testOrdering: \n" << *testOrdering; } else { cout << "Old inner product and new IP agree!!\n"; } } Teuchos::RCP<RHSEasy> rhs = Teuchos::rcp( new RHSEasy ); // the RHS as implemented by BurgersProblem divides the first component of beta by 2.0 // so we follow that. I've not done the math; just imitating the code... Teuchos::RCP<RHSEasy> otherRHS = Teuchos::rcp( new RHSEasy ); vector<double> e1_div2 = e1; e1_div2[0] /= 2.0; FunctionPtr rhsBeta = (e1_div2 * beta * e1 + Teuchos::rcp( new ConstantVectorFunction( e2 ) )) * u_prev; otherRHS->addTerm( rhsBeta * v->grad() - u_prev * tau->div() ); Teuchos::RCP<BurgersProblem> problem = Teuchos::rcp( new BurgersProblem(oldBurgersBF) ); expectedValues.resize(numCells, testOrdering->totalDofs() ); actualValues.resize (numCells, testOrdering->totalDofs() ); problem->integrateAgainstStandardBasis(expectedValues,testOrdering,basisCache); otherRHS->integrateAgainstStandardBasis(actualValues,testOrdering,basisCache); tol = 1e-14; maxDiff = 0.0; if (rank==0) { if ( ! TestSuite::fcsAgree(expectedValues,actualValues,tol,maxDiff) ) { cout << "Test failed: old RHS differs from new (\"otherRHS\"); maxDiff " << maxDiff << ".\n"; cout << "Old: \n" << expectedValues; cout << "New: \n" << actualValues; cout << "testOrdering: \n" << *testOrdering; } else { cout << "Old and new RHS (\"otherRHS\") agree!!\n"; } } FunctionPtr u_prev_squared_div2 = 0.5 * u_prev * u_prev; rhs->addTerm( (e1 * u_prev_squared_div2 + e2 * u_prev) * v->grad() - u_prev * tau->div()); if (! functionsAgree(e2 * u_prev, Teuchos::rcp( new ConstantVectorFunction( e2 ) ) * u_prev, basisCache) ) { cout << "two like functions differ...\n"; } FunctionPtr e1_f = Teuchos::rcp( new ConstantVectorFunction( e1 ) ); FunctionPtr e2_f = Teuchos::rcp( new ConstantVectorFunction( e2 ) ); FunctionPtr one = Teuchos::rcp( new ConstantScalarFunction( 1.0 ) ); if (! functionsAgree( Teuchos::rcp( new ProductFunction(e1_f, (e1_f + e2_f)) ), // e1_f * (e1_f + e2_f) one, basisCache) ) { cout << "two like functions differ...\n"; } if (! functionsAgree(u_prev_squared_div2, (e1_div2 * beta) * u_prev, basisCache) ) { cout << "two like functions differ...\n"; } if (! functionsAgree(e1 * u_prev_squared_div2, (e1_div2 * beta * e1) * u_prev, basisCache) ) { cout << "two like functions differ...\n"; } if (! functionsAgree(e1 * u_prev_squared_div2 + e2 * u_prev, (e1_div2 * beta * e1 + Teuchos::rcp( new ConstantVectorFunction( e2 ) )) * u_prev, basisCache) ) { cout << "two like functions differ...\n"; } problem->integrateAgainstStandardBasis(expectedValues,testOrdering,basisCache); rhs->integrateAgainstStandardBasis(actualValues,testOrdering,basisCache); tol = 1e-14; maxDiff = 0.0; if (rank==0) { if ( ! TestSuite::fcsAgree(expectedValues,actualValues,tol,maxDiff) ) { cout << "Test failed: old RHS differs from new (\"rhs\"); maxDiff " << maxDiff << ".\n"; cout << "Old: \n" << expectedValues; cout << "New: \n" << actualValues; cout << "testOrdering: \n" << *testOrdering; } else { cout << "Old and new RHS (\"rhs\") agree!!\n"; } } SpatialFilterPtr outflowBoundary = Teuchos::rcp( new TopBoundary ); SpatialFilterPtr inflowBoundary = Teuchos::rcp( new NegatedSpatialFilter(outflowBoundary) ); Teuchos::RCP<PenaltyConstraints> pc = Teuchos::rcp(new PenaltyConstraints); LinearTermPtr sigma_hat = beta * uhat->times_normal() - beta_n_u_minus_sigma_hat; FunctionPtr zero = Teuchos::rcp( new ConstantScalarFunction(0.0) ); pc->addConstraint(sigma_hat==zero,outflowBoundary); FunctionPtr u0 = Teuchos::rcp( new U0 ); FunctionPtr n = Teuchos::rcp( new UnitNormalFunction ); Teuchos::RCP<BCEasy> inflowBC = Teuchos::rcp( new BCEasy ); FunctionPtr u0_squared_div_2 = 0.5 * u0 * u0; inflowBC->addDirichlet(beta_n_u_minus_sigma_hat,inflowBoundary, ( e1 * u0_squared_div_2 + e2 * u0) * n ); // create a solution object Teuchos::RCP<Solution> solution = Teuchos::rcp(new Solution(mesh, inflowBC, rhs, ip)); mesh->registerSolution(solution); solution->setFilter(pc); // old penalty filter: Teuchos::RCP<LocalStiffnessMatrixFilter> penaltyBC = Teuchos::rcp(new PenaltyMethodFilter(problem)); // solution->setFilter(penaltyBC); // compare old and new filters elemType = mesh->getElement(0)->elementType(); trialOrdering = elemType->trialOrderPtr; testOrdering = elemType->testOrderPtr; // stiffness expectedValues.resize(numCells, trialOrdering->totalDofs(), trialOrdering->totalDofs() ); actualValues.resize (numCells, trialOrdering->totalDofs(), trialOrdering->totalDofs() ); expectedValues.initialize(0.0); actualValues.initialize(0.0); // load FieldContainer<double> expectedLoad(numCells, trialOrdering->totalDofs() ); FieldContainer<double> actualLoad(numCells, trialOrdering->totalDofs() ); penaltyBC->filter(expectedValues,expectedLoad,basisCache,mesh,problem); pc->filter(actualValues,actualLoad,basisCache,mesh,problem); maxDiff = 0.0; if (rank==0) { if ( ! TestSuite::fcsAgree(expectedValues,actualValues,tol,maxDiff) ) { cout << "Test failed: old penalty stiffness differs from new; maxDiff " << maxDiff << ".\n"; cout << "Old: \n" << expectedValues; cout << "New: \n" << actualValues; cout << "trialOrdering: \n" << *trialOrdering; } else { cout << "Old and new penalty stiffness agree!!\n"; } } if (rank==0) { if ( ! TestSuite::fcsAgree(expectedLoad,actualLoad,tol,maxDiff) ) { cout << "Test failed: old penalty load differs from new; maxDiff " << maxDiff << ".\n"; cout << "Old: \n" << expectedValues; cout << "New: \n" << actualValues; cout << "trialOrdering: \n" << *trialOrdering; } else { cout << "Old and new penalty load agree!!\n"; } } // define refinement strategy: Teuchos::RCP<RefinementStrategy> refinementStrategy = Teuchos::rcp(new RefinementStrategy(solution,energyThreshold)); // =================== END INITIALIZATION CODE ========================== // refine the spectral mesh, for comparability with the original Burgers' driver mesh->hRefine(vector<int>(1),RefinementPattern::regularRefinementPatternQuad()); int numRefs = 5; Teuchos::RCP<NonlinearStepSize> stepSize = Teuchos::rcp(new NonlinearStepSize(nonlinearStepSize)); Teuchos::RCP<NonlinearSolveStrategy> solveStrategy = Teuchos::rcp( new NonlinearSolveStrategy(backgroundFlow, solution, stepSize, nonlinearRelativeEnergyTolerance) ); for (int refIndex=0; refIndex<numRefs; refIndex++) { solveStrategy->solve(rank==0); refinementStrategy->refine(rank==0); // print to console on rank 0 } // one more nonlinear solve on refined mesh int numNRSteps = 5; for (int i=0; i<numNRSteps; i++) { solution->solve(false); backgroundFlow->addSolution(solution,1.0); } if (rank==0) { backgroundFlow->writeFieldsToFile(BurgersBilinearForm::U, "u_ref.m"); backgroundFlow->writeFieldsToFile(BurgersBilinearForm::SIGMA_1, "sigmax.m"); backgroundFlow->writeFieldsToFile(BurgersBilinearForm::SIGMA_2, "sigmay.m"); solution->writeFluxesToFile(BurgersBilinearForm::U_HAT, "du_hat_ref.dat"); } return 0; }
bool MeshTestUtility::neighborBasesAgreeOnSides(Teuchos::RCP<Mesh> mesh, Epetra_MultiVector &globalSolutionCoefficients) { bool success = true; MeshTopologyViewPtr meshTopo = mesh->getTopology(); int spaceDim = meshTopo->getDimension(); set<IndexType> activeCellIndices = meshTopo->getActiveCellIndices(); for (set<IndexType>::iterator cellIt=activeCellIndices.begin(); cellIt != activeCellIndices.end(); cellIt++) { IndexType cellIndex = *cellIt; CellPtr cell = meshTopo->getCell(cellIndex); BasisCachePtr fineCellBasisCache = BasisCache::basisCacheForCell(mesh, cellIndex); ElementTypePtr fineElemType = mesh->getElementType(cellIndex); DofOrderingPtr fineElemTrialOrder = fineElemType->trialOrderPtr; FieldContainer<double> fineSolutionCoefficients(fineElemTrialOrder->totalDofs()); mesh->globalDofAssignment()->interpretGlobalCoefficients(cellIndex, fineSolutionCoefficients, globalSolutionCoefficients); // if ((cellIndex==0) || (cellIndex==2)) { // cout << "MeshTestUtility: local coefficients for cell " << cellIndex << ":\n" << fineSolutionCoefficients; // } unsigned sideCount = cell->getSideCount(); for (unsigned sideOrdinal=0; sideOrdinal<sideCount; sideOrdinal++) { FieldContainer<double> fineSideRefPoints, fineCellRefPoints, coarseSideRefPoints, coarseCellRefPoints; bool hasCoarserNeighbor = determineRefTestPointsForNeighbors(meshTopo, cell, sideOrdinal, fineSideRefPoints, fineCellRefPoints, coarseSideRefPoints, coarseCellRefPoints); if (!hasCoarserNeighbor) continue; pair<GlobalIndexType, unsigned> neighborInfo = cell->getNeighborInfo(sideOrdinal, meshTopo); CellPtr neighborCell = meshTopo->getCell(neighborInfo.first); unsigned numTestPoints = coarseCellRefPoints.dimension(0); // cout << "testing neighbor agreement between cell " << cellIndex << " and " << neighborCell->cellIndex() << endl; // if we get here, the cell has a neighbor on this side, and is at least as fine as that neighbor. BasisCachePtr fineSideBasisCache = fineCellBasisCache->getSideBasisCache(sideOrdinal); fineCellBasisCache->setRefCellPoints(fineCellRefPoints); BasisCachePtr coarseCellBasisCache = BasisCache::basisCacheForCell(mesh, neighborInfo.first); BasisCachePtr coarseSideBasisCache = coarseCellBasisCache->getSideBasisCache(neighborInfo.second); coarseCellBasisCache->setRefCellPoints(coarseCellRefPoints); bool hasSidePoints = (fineSideRefPoints.size() > 0); if (hasSidePoints) { fineSideBasisCache->setRefCellPoints(fineSideRefPoints); coarseSideBasisCache->setRefCellPoints(coarseSideRefPoints); } // sanity check: do physical points match? FieldContainer<double> fineCellPhysicalCubaturePoints = fineCellBasisCache->getPhysicalCubaturePoints(); FieldContainer<double> coarseCellPhysicalCubaturePoints = coarseCellBasisCache->getPhysicalCubaturePoints(); double tol = 1e-14; double maxDiff = 0; if (! fcsAgree(coarseCellPhysicalCubaturePoints, fineCellPhysicalCubaturePoints, tol, maxDiff) ) { cout << "ERROR: MeshTestUtility::neighborBasesAgreeOnSides internal error: fine and coarse cell cubature points do not agree.\n"; success = false; continue; } if (hasSidePoints) { FieldContainer<double> fineSidePhysicalCubaturePoints = fineSideBasisCache->getPhysicalCubaturePoints(); FieldContainer<double> coarseSidePhysicalCubaturePoints = coarseSideBasisCache->getPhysicalCubaturePoints(); double tol = 1e-14; double maxDiff = 0; if (! fcsAgree(fineSidePhysicalCubaturePoints, fineCellPhysicalCubaturePoints, tol, maxDiff) ) { cout << "ERROR: MeshTestUtility::neighborBasesAgreeOnSides internal error: fine side and cell cubature points do not agree.\n"; success = false; continue; } if (! fcsAgree(coarseSidePhysicalCubaturePoints, coarseCellPhysicalCubaturePoints, tol, maxDiff) ) { cout << "ERROR: MeshTestUtility::neighborBasesAgreeOnSides internal error: coarse side and cell cubature points do not agree.\n"; success = false; continue; } if (! fcsAgree(coarseSidePhysicalCubaturePoints, fineSidePhysicalCubaturePoints, tol, maxDiff) ) { cout << "ERROR: MeshTestUtility::neighborBasesAgreeOnSides internal error: fine and coarse side cubature points do not agree.\n"; success = false; continue; } } ElementTypePtr coarseElementType = mesh->getElementType(neighborInfo.first); DofOrderingPtr coarseElemTrialOrder = coarseElementType->trialOrderPtr; FieldContainer<double> coarseSolutionCoefficients(coarseElemTrialOrder->totalDofs()); mesh->globalDofAssignment()->interpretGlobalCoefficients(neighborInfo.first, coarseSolutionCoefficients, globalSolutionCoefficients); set<int> varIDs = fineElemTrialOrder->getVarIDs(); for (set<int>::iterator varIt = varIDs.begin(); varIt != varIDs.end(); varIt++) { int varID = *varIt; // cout << "MeshTestUtility: varID " << varID << ":\n"; bool isTraceVar = mesh->bilinearForm()->isFluxOrTrace(varID); BasisPtr fineBasis, coarseBasis; BasisCachePtr fineCache, coarseCache; if (isTraceVar) { if (! hasSidePoints) continue; // then nothing to do for traces fineBasis = fineElemTrialOrder->getBasis(varID, sideOrdinal); coarseBasis = coarseElemTrialOrder->getBasis(varID, neighborInfo.second); fineCache = fineSideBasisCache; coarseCache = coarseSideBasisCache; } else { fineBasis = fineElemTrialOrder->getBasis(varID); coarseBasis = coarseElemTrialOrder->getBasis(varID); fineCache = fineCellBasisCache; coarseCache = coarseCellBasisCache; Camellia::EFunctionSpace fs = fineBasis->functionSpace(); if ((fs == Camellia::FUNCTION_SPACE_HVOL) || (fs == Camellia::FUNCTION_SPACE_VECTOR_HVOL) || (fs == Camellia::FUNCTION_SPACE_TENSOR_HVOL)) { // volume L^2 basis: no continuities expected... continue; } } FieldContainer<double> localValuesFine = *fineCache->getTransformedValues(fineBasis, OP_VALUE); FieldContainer<double> localValuesCoarse = *coarseCache->getTransformedValues(coarseBasis, OP_VALUE); bool scalarValued = (localValuesFine.rank() == 3); // the following used if vector or tensor-valued: Teuchos::Array<int> valueDim; unsigned componentsPerValue = 1; FieldContainer<double> valueContainer; // just used for enumeration computation... if (!scalarValued) { localValuesFine.dimensions(valueDim); // clear first three: valueDim.erase(valueDim.begin()); valueDim.erase(valueDim.begin()); valueDim.erase(valueDim.begin()); valueContainer.resize(valueDim); componentsPerValue = valueContainer.size(); } // if (localValuesFine.rank() != 3) { // cout << "WARNING: MeshTestUtility::neighborBasesAgreeOnSides() only supports scalar-valued bases right now. Skipping check for varID " << varID << endl; // continue; // } FieldContainer<double> localPointValuesFine(fineElemTrialOrder->totalDofs()); FieldContainer<double> localPointValuesCoarse(coarseElemTrialOrder->totalDofs()); for (int valueComponentOrdinal=0; valueComponentOrdinal<componentsPerValue; valueComponentOrdinal++) { Teuchos::Array<int> valueMultiIndex(valueContainer.rank()); if (!scalarValued) valueContainer.getMultiIndex(valueMultiIndex, valueComponentOrdinal); Teuchos::Array<int> localValuesMultiIndex(localValuesFine.rank()); for (int r=0; r<valueMultiIndex.size(); r++) { localValuesMultiIndex[r+3] = valueMultiIndex[r]; } for (int ptOrdinal=0; ptOrdinal<numTestPoints; ptOrdinal++) { localPointValuesCoarse.initialize(0); localPointValuesFine.initialize(0); localValuesMultiIndex[2] = ptOrdinal; double fineSolutionValue = 0, coarseSolutionValue = 0; for (int basisOrdinal=0; basisOrdinal < fineBasis->getCardinality(); basisOrdinal++) { int fineDofIndex; if (isTraceVar) fineDofIndex = fineElemTrialOrder->getDofIndex(varID, basisOrdinal, sideOrdinal); else fineDofIndex = fineElemTrialOrder->getDofIndex(varID, basisOrdinal); if (scalarValued) { localPointValuesFine(fineDofIndex) = localValuesFine(0,basisOrdinal,ptOrdinal); } else { localValuesMultiIndex[1] = basisOrdinal; localPointValuesFine(fineDofIndex) = localValuesFine.getValue(localValuesMultiIndex); } fineSolutionValue += fineSolutionCoefficients(fineDofIndex) * localPointValuesFine(fineDofIndex); } for (int basisOrdinal=0; basisOrdinal < coarseBasis->getCardinality(); basisOrdinal++) { int coarseDofIndex; if (isTraceVar) coarseDofIndex = coarseElemTrialOrder->getDofIndex(varID, basisOrdinal, neighborInfo.second); else coarseDofIndex = coarseElemTrialOrder->getDofIndex(varID, basisOrdinal); if (scalarValued) { localPointValuesCoarse(coarseDofIndex) = localValuesCoarse(0,basisOrdinal,ptOrdinal); } else { localValuesMultiIndex[1] = basisOrdinal; localPointValuesCoarse(coarseDofIndex) = localValuesCoarse.getValue(localValuesMultiIndex); } coarseSolutionValue += coarseSolutionCoefficients(coarseDofIndex) * localPointValuesCoarse(coarseDofIndex); } if (abs(coarseSolutionValue - fineSolutionValue) > 1e-13) { success = false; cout << "coarseSolutionValue (" << coarseSolutionValue << ") and fineSolutionValue (" << fineSolutionValue << ") differ by " << abs(coarseSolutionValue - fineSolutionValue); cout << " at point " << ptOrdinal << " for varID " << varID << ". "; cout << "This may be an indication that something is amiss with the global-to-local map.\n"; } else { // // DEBUGGING: // cout << "solution value at point ("; // for (int d=0; d<spaceDim-1; d++) { // cout << fineSidePhysicalCubaturePoints(0,ptOrdinal,d) << ", "; // } // cout << fineSidePhysicalCubaturePoints(0,ptOrdinal,spaceDim-1) << "): "; // cout << coarseSolutionValue << endl; } FieldContainer<double> globalValuesFromFine, globalValuesFromCoarse; FieldContainer<GlobalIndexType> globalDofIndicesFromFine, globalDofIndicesFromCoarse; mesh->globalDofAssignment()->interpretLocalData(cellIndex, localPointValuesFine, globalValuesFromFine, globalDofIndicesFromFine); mesh->globalDofAssignment()->interpretLocalData(neighborInfo.first, localPointValuesCoarse, globalValuesFromCoarse, globalDofIndicesFromCoarse); std::map<GlobalIndexType, double> fineValuesMap; std::map<GlobalIndexType, double> coarseValuesMap; for (int i=0; i<globalDofIndicesFromCoarse.size(); i++) { GlobalIndexType globalDofIndex = globalDofIndicesFromCoarse[i]; coarseValuesMap[globalDofIndex] = globalValuesFromCoarse[i]; } double maxDiff = 0; for (int i=0; i<globalDofIndicesFromFine.size(); i++) { GlobalIndexType globalDofIndex = globalDofIndicesFromFine[i]; fineValuesMap[globalDofIndex] = globalValuesFromFine[i]; double diff = abs( fineValuesMap[globalDofIndex] - coarseValuesMap[globalDofIndex]); maxDiff = std::max(diff, maxDiff); if (diff > tol) { success = false; cout << "interpreted fine and coarse disagree at point ("; for (int d=0; d<spaceDim; d++) { cout << fineCellPhysicalCubaturePoints(0,ptOrdinal,d); if (d==spaceDim-1) cout << ").\n"; else cout << ", "; } } } if (maxDiff > tol) { cout << "maxDiff: " << maxDiff << endl; cout << "globalValuesFromFine:\n" << globalValuesFromFine; cout << "globalValuesFromCoarse:\n" << globalValuesFromCoarse; cout << "globalDofIndicesFromFine:\n" << globalDofIndicesFromFine; cout << "globalDofIndicesFromCoarse:\n" << globalDofIndicesFromCoarse; continue; // only worth testing further if we passed the above } } } } } } // cout << "Completed neighborBasesAgreeOnSides.\n"; return success; }