bool ScratchPadTests::testConstantFunctionProduct() { bool success = true; // set up basisCache (even though it won't really be used here) BasisCachePtr basisCache = Teuchos::rcp( new BasisCache( _elemType, _spectralConfusionMesh ) ); vector<GlobalIndexType> cellIDs; int cellID = 0; cellIDs.push_back(cellID); basisCache->setPhysicalCellNodes( _spectralConfusionMesh->physicalCellNodesForCell(cellID), cellIDs, true ); int numCells = _basisCache->getPhysicalCubaturePoints().dimension(0); int numPoints = _testPoints.dimension(0); FunctionPtr three = Function::constant(3.0); FunctionPtr two = Function::constant(2.0); FieldContainer<double> values(numCells,numPoints); two->values(values,basisCache); three->scalarMultiplyBasisValues( values, basisCache ); FieldContainer<double> expectedValues(numCells,numPoints); expectedValues.initialize( 3.0 * 2.0 ); double tol = 1e-15; double maxDiff = 0.0; if ( ! fcsAgree(expectedValues, values, tol, maxDiff) ) { success = false; cout << "Expected product differs from actual; maxDiff: " << maxDiff << endl; } return success; }
bool FunctionTests::testAdaptiveIntegrate() { bool success = true; // we must create our own basisCache here because _basisCache // has had its ref cell points set, which basically means it's // opted out of having any help with integration. BasisCachePtr basisCache = Teuchos::rcp( new BasisCache( _elemType, _spectralConfusionMesh ) ); vector<GlobalIndexType> cellIDs; cellIDs.push_back(0); basisCache->setPhysicalCellNodes( _spectralConfusionMesh->physicalCellNodesForCell(0), cellIDs, true ); double eps = .1; // FunctionPtr boundaryLayerFunction = Teuchos::rcp( new BoundaryLayerFunction(eps) ); int numCells = basisCache->cellIDs().size(); FieldContainer<double> integrals(numCells); double quadtol = 1e-2; double computedIntegral = boundaryLayerFunction->integrate(_spectralConfusionMesh,quadtol); double trueIntegral = (eps*(exp(1/eps) - exp(-1/eps)))*2.0; double diff = trueIntegral-computedIntegral; double relativeError = abs(diff)/abs(trueIntegral); // relative error double tol = 1e-2; if (relativeError > tol) { success = false; cout << "failing testAdaptiveIntegrate() with computed integral " << computedIntegral << " and true integral " << trueIntegral << endl; } return success; }
void PreviousSolutionFunction::values(FieldContainer<double> &values, BasisCachePtr basisCache) { int rank = Teuchos::GlobalMPISession::getRank(); if (_overrideMeshCheck) { _solnExpression->evaluate(values, _soln, basisCache); return; } if (!basisCache.get()) cout << "basisCache is nil!\n"; if (!_soln.get()) cout << "_soln is nil!\n"; // values are stored in (C,P,D) order if (basisCache->mesh().get() == _soln->mesh().get()) { _solnExpression->evaluate(values, _soln, basisCache); } else { static bool warningIssued = false; if (!warningIssued) { if (rank==0) cout << "NOTE: In PreviousSolutionFunction, basisCache's mesh doesn't match solution's. If this is not what you intended, it would be a good idea to make sure that the mesh is passed in on BasisCache construction; the evaluation will be a lot slower without it...\n"; warningIssued = true; } // get the physicalPoints, and make a basisCache for each... FieldContainer<double> physicalPoints = basisCache->getPhysicalCubaturePoints(); FieldContainer<double> value(1,1); // assumes scalar-valued solution function. int numCells = physicalPoints.dimension(0); int numPoints = physicalPoints.dimension(1); int spaceDim = physicalPoints.dimension(2); physicalPoints.resize(numCells*numPoints,spaceDim); vector< ElementPtr > elements = _soln->mesh()->elementsForPoints(physicalPoints, false); // false: don't make elements null just because they're off-rank. FieldContainer<double> point(1,spaceDim); FieldContainer<double> refPoint(1,spaceDim); int combinedIndex = 0; vector<GlobalIndexType> cellID; cellID.push_back(-1); BasisCachePtr basisCacheOnePoint; for (int cellIndex=0; cellIndex<numCells; cellIndex++) { for (int ptIndex=0; ptIndex<numPoints; ptIndex++, combinedIndex++) { if (elements[combinedIndex].get()==NULL) continue; // no element found for point; skip it… ElementTypePtr elemType = elements[combinedIndex]->elementType(); for (int d=0; d<spaceDim; d++) { point(0,d) = physicalPoints(combinedIndex,d); } if (elements[combinedIndex]->cellID() != cellID[0]) { cellID[0] = elements[combinedIndex]->cellID(); basisCacheOnePoint = Teuchos::rcp( new BasisCache(elemType, _soln->mesh()) ); basisCacheOnePoint->setPhysicalCellNodes(_soln->mesh()->physicalCellNodesForCell(cellID[0]),cellID,false); // false: don't createSideCacheToo } // compute the refPoint: typedef CellTools<double> CellTools; int whichCell = 0; CellTools::mapToReferenceFrame(refPoint,point,_soln->mesh()->physicalCellNodesForCell(cellID[0]), *(elemType->cellTopoPtr),whichCell); basisCacheOnePoint->setRefCellPoints(refPoint); // cout << "refCellPoints:\n " << refPoint; // cout << "physicalCubaturePoints:\n " << basisCacheOnePoint->getPhysicalCubaturePoints(); _solnExpression->evaluate(value, _soln, basisCacheOnePoint); // cout << "value at point (" << point(0,0) << ", " << point(0,1) << ") = " << value(0,0) << endl; values(cellIndex,ptIndex) = value(0,0); } } } }
bool RHSTests::testIntegrateAgainstStandardBasis() { bool success = true; double tol = 1e-14; Teuchos::RCP<ElementType> elemType = _mesh->getElement(0)->elementType(); int rank = _mesh->Comm()->MyPID(); vector< Teuchos::RCP< Element > > elemsInPartitionOfType = _mesh->elementsOfType(rank, elemType); FieldContainer<double> physicalCellNodes = _mesh->physicalCellNodes(elemType); int numCells = elemsInPartitionOfType.size(); int numTestDofs = elemType->testOrderPtr->totalDofs(); // set up diagonal testWeights matrices so we can reuse the existing computeRHS, and compare results… FieldContainer<double> testWeights(numCells,numTestDofs,numTestDofs); for (int cellIndex=0; cellIndex<numCells; cellIndex++) { for (int i=0; i<numTestDofs; i++) { testWeights(cellIndex,i,i) = 1.0; } } FieldContainer<double> rhsExpected(numCells,numTestDofs); FieldContainer<double> rhsActual(numCells,numTestDofs); // determine cellIDs vector<GlobalIndexType> cellIDs = _mesh->globalDofAssignment()->cellIDsOfElementType(rank, elemType); if (numCells > 0) { // prepare basisCache and cellIDs BasisCachePtr basisCache = Teuchos::rcp(new BasisCache(elemType,_mesh)); bool createSideCacheToo = true; basisCache->setPhysicalCellNodes(physicalCellNodes,cellIDs,createSideCacheToo); _rhs->integrateAgainstStandardBasis(rhsActual, elemType->testOrderPtr, basisCache); _rhs->integrateAgainstOptimalTests(rhsExpected, testWeights, elemType->testOrderPtr, basisCache); } double maxDiff = 0.0; if ( ! fcsAgree(rhsExpected,rhsActual,tol,maxDiff) ) { success = false; cout << "Failed testIntegrateAgainstStandardBasis: maxDiff = " << maxDiff << endl; } // check success across MPI nodes return allSuccess(success); }
bool RHSTests::testTrivialRHS() { bool success = true; double tol = 1e-14; int rank = _mesh->Comm()->MyPID(); Teuchos::RCP<ElementType> elemType = _mesh->getElement(0)->elementType(); vector< Teuchos::RCP< Element > > elemsInPartitionOfType = _mesh->elementsOfType(rank, elemType); FieldContainer<double> physicalCellNodes = _mesh->physicalCellNodes(elemType); int numCells = elemsInPartitionOfType.size(); int numTestDofs = elemType->testOrderPtr->totalDofs(); // determine cellIDs vector<GlobalIndexType> cellIDs = _mesh->globalDofAssignment()->cellIDsOfElementType(rank, elemType); if (numCells > 0) { // prepare basisCache and cellIDs BasisCachePtr basisCache = Teuchos::rcp(new BasisCache(elemType,_mesh)); bool createSideCacheToo = true; basisCache->setPhysicalCellNodes(physicalCellNodes,cellIDs,createSideCacheToo); FieldContainer<double> rhsExpected(numCells,numTestDofs); FunctionPtr zero = Function::constant(0.0); RHSPtr rhs = RHS::rhs(); FunctionPtr f = zero; rhs->addTerm( f * _v ); // obviously, with f = 0 adding this term is not necessary! rhs->integrateAgainstStandardBasis(rhsExpected, elemType->testOrderPtr, basisCache); for (int i = 0; i<numCells; i++) { for (int j = 0; j<numTestDofs; j++) { if (abs(rhsExpected(i,j))>tol) { success = false; } } } } return allSuccess(success); }
bool RHSTests::testRHSEasy() { bool success = true; double tol = 1e-14; int rank = _mesh->Comm()->MyPID(); int numProcs = Teuchos::GlobalMPISession::getNProc(); Teuchos::RCP<ElementType> elemType = _mesh->getElement(0)->elementType(); vector< Teuchos::RCP< Element > > elemsInPartitionOfType = _mesh->elementsOfType(rank, elemType); FieldContainer<double> physicalCellNodes = _mesh->physicalCellNodes(elemType); int numCells = elemsInPartitionOfType.size(); int numTestDofs = elemType->testOrderPtr->totalDofs(); FieldContainer<double> rhsExpected(numCells,numTestDofs); FieldContainer<double> rhsActual(numCells,numTestDofs); // determine cellIDs vector<GlobalIndexType> cellIDs = _mesh->globalDofAssignment()->cellIDsOfElementType(rank, elemType); // prepare basisCache and cellIDs if (numCells > 0) { BasisCachePtr basisCache = Teuchos::rcp(new BasisCache(elemType,_mesh)); bool createSideCacheToo = true; basisCache->setPhysicalCellNodes(physicalCellNodes,cellIDs,createSideCacheToo); _rhs->integrateAgainstStandardBasis(rhsActual, elemType->testOrderPtr, basisCache); _rhsEasy->integrateAgainstStandardBasis(rhsExpected, elemType->testOrderPtr, basisCache); } double maxDiff = 0.0; if ( ! fcsAgree(rhsExpected,rhsActual,tol,maxDiff) ) { success = false; cout << "Failed testRHSEasy: maxDiff = " << maxDiff << endl; cout << "Expected values:\n" << rhsExpected; cout << "Actual values:\n" << rhsActual; cout << "Test dof ordering:\n" << *(elemType->testOrderPtr); } return allSuccess(success); }
bool LinearTermTests::testIntegrateMixedBasis() { bool success = true; //////////////////// DECLARE VARIABLES /////////////////////// // define test variables VarFactoryPtr varFactory = VarFactory::varFactory(); VarPtr v = varFactory->testVar("v", HGRAD); // define trial variables VarPtr beta_n_u_hat = varFactory->fluxVar("\\widehat{\\beta \\cdot n }"); VarPtr u = varFactory->fieldVar("u"); vector<double> beta; beta.push_back(1.0); beta.push_back(1.0); //////////////////// DEFINE BILINEAR FORM/Mesh /////////////////////// BFPtr convectionBF = Teuchos::rcp( new BF(varFactory) ); // v terms: convectionBF->addTerm( -u, beta * v->grad() ); convectionBF->addTerm( beta_n_u_hat, v); convectionBF->addTerm( u, v); // build CONSTANT SINGLE ELEMENT MESH int order = 0; int H1Order = order+1; int pToAdd = 1; int nCells = 2; // along a side // create a pointer to a new mesh: Teuchos::RCP<Mesh> mesh = MeshUtilities::buildUnitQuadMesh(nCells,convectionBF, H1Order, H1Order+pToAdd); ElementTypePtr elemType = mesh->getElement(0)->elementType(); BasisCachePtr basisCache = Teuchos::rcp(new BasisCache(elemType, mesh)); vector<GlobalIndexType> cellIDs; vector< ElementPtr > allElems = mesh->activeElements(); vector< ElementPtr >::iterator elemIt; for (elemIt=allElems.begin(); elemIt!=allElems.end(); elemIt++) { cellIDs.push_back((*elemIt)->cellID()); } bool createSideCacheToo = true; basisCache->setPhysicalCellNodes(mesh->physicalCellNodesGlobal(elemType), cellIDs, createSideCacheToo); int numTrialDofs = elemType->trialOrderPtr->totalDofs(); int numCells = mesh->numActiveElements(); double areaPerCell = 1.0 / numCells; FieldContainer<double> integrals(numCells,numTrialDofs); FieldContainer<double> expectedIntegrals(numCells,numTrialDofs); double sidelengthOfCell = 1.0 / nCells; DofOrderingPtr trialOrdering = elemType->trialOrderPtr; int dofForField = trialOrdering->getDofIndex(u->ID(), 0); vector<int> dofsForFlux; const vector<int>* sidesForFlux = &trialOrdering->getSidesForVarID(beta_n_u_hat->ID()); for (vector<int>::const_iterator sideIt = sidesForFlux->begin(); sideIt != sidesForFlux->end(); sideIt++) { int sideIndex = *sideIt; dofsForFlux.push_back(trialOrdering->getDofIndex(beta_n_u_hat->ID(), 0, sideIndex)); } for (int cellIndex = 0; cellIndex < numCells; cellIndex++) { expectedIntegrals(cellIndex, dofForField) = areaPerCell; for (vector<int>::iterator dofIt = dofsForFlux.begin(); dofIt != dofsForFlux.end(); dofIt++) { int fluxDofIndex = *dofIt; expectedIntegrals(cellIndex, fluxDofIndex) = sidelengthOfCell; } } // cout << "expectedIntegrals:\n" << expectedIntegrals; // setup: with constant degrees of freedom, we expect that the integral of int_dK (flux) + int_K (field) will be ones for each degree of freedom, assuming the basis functions for these constants field/flux variables are just C = 1.0. // //On a unit square, int_K (constant) = 1.0, and int_dK (u_i) = 1, for i = 0,...,3. LinearTermPtr lt = 1.0 * beta_n_u_hat; LinearTermPtr field = 1.0 * u; lt->addTerm(field,true); lt->integrate(integrals, elemType->trialOrderPtr, basisCache); double tol = 1e-12; double maxDiff; success = TestSuite::fcsAgree(integrals,expectedIntegrals,tol,maxDiff); if (success==false) { cout << "Failed testIntegrateMixedBasis with maxDiff = " << maxDiff << endl; } return success; }
int main(int argc, char *argv[]) { int rank = 0; #ifdef HAVE_MPI // TODO: figure out the right thing to do here... // may want to modify argc and argv before we make the following call: Teuchos::GlobalMPISession mpiSession(&argc, &argv,0); rank=mpiSession.getRank(); #else #endif bool useLineSearch = false; int pToAdd = 2; // for optimal test function approximation int pToAddForStreamFunction = 2; double nonlinearStepSize = 1.0; double dt = 0.5; double nonlinearRelativeEnergyTolerance = 0.015; // used to determine convergence of the nonlinear solution // double nonlinearRelativeEnergyTolerance = 0.15; // used to determine convergence of the nonlinear solution double eps = 1.0/64.0; // width of ramp up to 1.0 for top BC; eps == 0 ==> soln not in H1 // epsilon above is chosen to match our initial 16x16 mesh, to avoid quadrature errors. // double eps = 0.0; // John Evans's problem: not in H^1 bool enforceLocalConservation = false; bool enforceOneIrregularity = true; bool reportPerCellErrors = true; bool useMumps = true; int horizontalCells, verticalCells; int maxIters = 50; // for nonlinear steps vector<double> ReValues; // usage: polyOrder [numRefinements] // parse args: if (argc < 6) { cout << "Usage: NavierStokesCavityFlowContinuationFixedMesh fieldPolyOrder hCells vCells energyErrorGoal Re0 [Re1 ...]\n"; return -1; } int polyOrder = atoi(argv[1]); horizontalCells = atoi(argv[2]); verticalCells = atoi(argv[3]); double energyErrorGoal = atof(argv[4]); for (int i=5; i<argc; i++) { ReValues.push_back(atof(argv[i])); } if (rank == 0) { cout << "L^2 order: " << polyOrder << endl; cout << "initial mesh size: " << horizontalCells << " x " << verticalCells << endl; cout << "energy error goal: " << energyErrorGoal << endl; cout << "Reynolds number values for continuation:\n"; for (int i=0; i<ReValues.size(); i++) { cout << ReValues[i] << ", "; } cout << endl; } 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; // define meshes: int H1Order = polyOrder + 1; bool useTriangles = false; bool meshHasTriangles = useTriangles; double minL2Increment = 1e-8; // get variable definitions: VarFactory varFactory = VGPStokesFormulation::vgpVarFactory(); u1 = varFactory.fieldVar(VGP_U1_S); u2 = varFactory.fieldVar(VGP_U2_S); sigma11 = varFactory.fieldVar(VGP_SIGMA11_S); sigma12 = varFactory.fieldVar(VGP_SIGMA12_S); sigma21 = varFactory.fieldVar(VGP_SIGMA21_S); sigma22 = varFactory.fieldVar(VGP_SIGMA22_S); p = varFactory.fieldVar(VGP_P_S); u1hat = varFactory.traceVar(VGP_U1HAT_S); u2hat = varFactory.traceVar(VGP_U2HAT_S); t1n = varFactory.fluxVar(VGP_T1HAT_S); t2n = varFactory.fluxVar(VGP_T2HAT_S); v1 = varFactory.testVar(VGP_V1_S, HGRAD); v2 = varFactory.testVar(VGP_V2_S, HGRAD); tau1 = varFactory.testVar(VGP_TAU1_S, HDIV); tau2 = varFactory.testVar(VGP_TAU2_S, HDIV); q = varFactory.testVar(VGP_Q_S, HGRAD); FunctionPtr u1_0 = Teuchos::rcp( new U1_0(eps) ); FunctionPtr u2_0 = Teuchos::rcp( new U2_0 ); FunctionPtr zero = Function::zero(); ParameterFunctionPtr Re_param = ParameterFunction::parameterFunction(1); VGPNavierStokesProblem problem = VGPNavierStokesProblem(Re_param,quadPoints, horizontalCells,verticalCells, H1Order, pToAdd, u1_0, u2_0, // BC for u zero, zero); // zero forcing function SolutionPtr solution = problem.backgroundFlow(); SolutionPtr solnIncrement = problem.solutionIncrement(); Teuchos::RCP<Mesh> mesh = problem.mesh(); mesh->registerSolution(solution); mesh->registerSolution(solnIncrement); /////////////////////////////////////////////////////////////////////////// // define bilinear form for stream function: VarFactory streamVarFactory; VarPtr phi_hat = streamVarFactory.traceVar("\\widehat{\\phi}"); VarPtr psin_hat = streamVarFactory.fluxVar("\\widehat{\\psi}_n"); VarPtr psi_1 = streamVarFactory.fieldVar("\\psi_1"); VarPtr psi_2 = streamVarFactory.fieldVar("\\psi_2"); VarPtr phi = streamVarFactory.fieldVar("\\phi"); VarPtr q_s = streamVarFactory.testVar("q_s", HGRAD); VarPtr v_s = streamVarFactory.testVar("v_s", HDIV); BFPtr streamBF = Teuchos::rcp( new BF(streamVarFactory) ); streamBF->addTerm(psi_1, q_s->dx()); streamBF->addTerm(psi_2, q_s->dy()); streamBF->addTerm(-psin_hat, q_s); streamBF->addTerm(psi_1, v_s->x()); streamBF->addTerm(psi_2, v_s->y()); streamBF->addTerm(phi, v_s->div()); streamBF->addTerm(-phi_hat, v_s->dot_normal()); Teuchos::RCP<Mesh> streamMesh, overkillMesh; streamMesh = MeshFactory::buildQuadMesh(quadPoints, horizontalCells, verticalCells, streamBF, H1Order+pToAddForStreamFunction, H1Order+pToAdd+pToAddForStreamFunction, useTriangles); mesh->registerObserver(streamMesh); // will refine streamMesh in the same way as mesh. map<int, double> dofsToL2error; // key: numGlobalDofs, value: total L2error compared with overkill vector< VarPtr > fields; fields.push_back(u1); fields.push_back(u2); fields.push_back(sigma11); fields.push_back(sigma12); fields.push_back(sigma21); fields.push_back(sigma22); fields.push_back(p); if (rank == 0) { cout << "Starting mesh has " << horizontalCells << " x " << verticalCells << " elements and "; cout << mesh->numGlobalDofs() << " total dofs.\n"; cout << "polyOrder = " << polyOrder << endl; cout << "pToAdd = " << pToAdd << endl; cout << "eps for top BC = " << eps << endl; if (useTriangles) { cout << "Using triangles.\n"; } if (enforceLocalConservation) { cout << "Enforcing local conservation.\n"; } else { cout << "NOT enforcing local conservation.\n"; } if (enforceOneIrregularity) { cout << "Enforcing 1-irregularity.\n"; } else { cout << "NOT enforcing 1-irregularity.\n"; } } //////////////////// CREATE BCs /////////////////////// SpatialFilterPtr entireBoundary = Teuchos::rcp( new SpatialFilterUnfiltered ); FunctionPtr u1_prev = Function::solution(u1,solution); FunctionPtr u2_prev = Function::solution(u2,solution); FunctionPtr u1hat_prev = Function::solution(u1hat,solution); FunctionPtr u2hat_prev = Function::solution(u2hat,solution); //////////////////// SOLVE & REFINE /////////////////////// FunctionPtr vorticity = Teuchos::rcp( new PreviousSolutionFunction(solution, - u1->dy() + u2->dx() ) ); // FunctionPtr vorticity = Teuchos::rcp( new PreviousSolutionFunction(solution,sigma12 - sigma21) ); RHSPtr streamRHS = RHS::rhs(); streamRHS->addTerm(vorticity * q_s); ((PreviousSolutionFunction*) vorticity.get())->setOverrideMeshCheck(true); ((PreviousSolutionFunction*) u1_prev.get())->setOverrideMeshCheck(true); ((PreviousSolutionFunction*) u2_prev.get())->setOverrideMeshCheck(true); BCPtr streamBC = BC::bc(); // streamBC->addDirichlet(psin_hat, entireBoundary, u0_cross_n); streamBC->addDirichlet(phi_hat, entireBoundary, zero); // streamBC->addZeroMeanConstraint(phi); IPPtr streamIP = Teuchos::rcp( new IP ); streamIP->addTerm(q_s); streamIP->addTerm(q_s->grad()); streamIP->addTerm(v_s); streamIP->addTerm(v_s->div()); SolutionPtr streamSolution = Teuchos::rcp( new Solution( streamMesh, streamBC, streamRHS, streamIP ) ); if (enforceLocalConservation) { FunctionPtr zero = Function::zero(); solution->lagrangeConstraints()->addConstraint(u1hat->times_normal_x() + u2hat->times_normal_y()==zero); solnIncrement->lagrangeConstraints()->addConstraint(u1hat->times_normal_x() + u2hat->times_normal_y()==zero); } if (true) { FunctionPtr u1_incr = Function::solution(u1, solnIncrement); FunctionPtr u2_incr = Function::solution(u2, solnIncrement); FunctionPtr sigma11_incr = Function::solution(sigma11, solnIncrement); FunctionPtr sigma12_incr = Function::solution(sigma12, solnIncrement); FunctionPtr sigma21_incr = Function::solution(sigma21, solnIncrement); FunctionPtr sigma22_incr = Function::solution(sigma22, solnIncrement); FunctionPtr p_incr = Function::solution(p, solnIncrement); FunctionPtr l2_incr = u1_incr * u1_incr + u2_incr * u2_incr + p_incr * p_incr + sigma11_incr * sigma11_incr + sigma12_incr * sigma12_incr + sigma21_incr * sigma21_incr + sigma22_incr * sigma22_incr; double energyThreshold = 0.20; Teuchos::RCP< RefinementStrategy > refinementStrategy = Teuchos::rcp( new RefinementStrategy( solnIncrement, energyThreshold )); for (int i=0; i<ReValues.size(); i++) { double Re = ReValues[i]; Re_param->setValue(Re); if (rank==0) cout << "Solving with Re = " << Re << ":\n"; double energyErrorTotal; do { double incr_norm; do { problem.iterate(useLineSearch); incr_norm = sqrt(l2_incr->integrate(problem.mesh())); if (rank==0) { cout << "\x1B[2K"; // Erase the entire current line. cout << "\x1B[0E"; // Move to the beginning of the current line. cout << "Iteration: " << problem.iterationCount() << "; L^2(incr) = " << incr_norm; flush(cout); } } while ((incr_norm > minL2Increment ) && (problem.iterationCount() < maxIters)); if (rank==0) cout << endl; problem.setIterationCount(1); // 1 means reuse background flow (which we must, given that we want continuation in Re...) energyErrorTotal = solnIncrement->energyErrorTotal(); //solution->energyErrorTotal(); if (energyErrorTotal > energyErrorGoal) { refinementStrategy->refine(false); } if (rank==0) { cout << "Energy error: " << energyErrorTotal << endl; } } while (energyErrorTotal > energyErrorGoal); } } double energyErrorTotal = solution->energyErrorTotal(); double incrementalEnergyErrorTotal = solnIncrement->energyErrorTotal(); if (rank == 0) { cout << "final mesh has " << mesh->numActiveElements() << " elements and " << mesh->numGlobalDofs() << " dofs.\n"; cout << "energy error: " << energyErrorTotal << endl; cout << " (Incremental solution's energy error is " << incrementalEnergyErrorTotal << ".)\n"; } FunctionPtr u1_sq = u1_prev * u1_prev; FunctionPtr u_dot_u = u1_sq + (u2_prev * u2_prev); FunctionPtr u_mag = Teuchos::rcp( new SqrtFunction( u_dot_u ) ); FunctionPtr u_div = Teuchos::rcp( new PreviousSolutionFunction(solution, u1->dx() + u2->dy() ) ); FunctionPtr massFlux = Teuchos::rcp( new PreviousSolutionFunction(solution, u1hat->times_normal_x() + u2hat->times_normal_y()) ); // check that the zero mean pressure is being correctly imposed: FunctionPtr p_prev = Teuchos::rcp( new PreviousSolutionFunction(solution,p) ); double p_avg = p_prev->integrate(mesh); if (rank==0) cout << "Integral of pressure: " << p_avg << endl; // integrate massFlux over each element (a test): // fake a new bilinear form so we can integrate against 1 VarPtr testOne = varFactory.testVar("1",CONSTANT_SCALAR); BFPtr fakeBF = Teuchos::rcp( new BF(varFactory) ); LinearTermPtr massFluxTerm = massFlux * testOne; CellTopoPtrLegacy quadTopoPtr = Teuchos::rcp(new shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() )); DofOrderingFactory dofOrderingFactory(fakeBF); int fakeTestOrder = H1Order; DofOrderingPtr testOrdering = dofOrderingFactory.testOrdering(fakeTestOrder, *quadTopoPtr); int testOneIndex = testOrdering->getDofIndex(testOne->ID(),0); vector< ElementTypePtr > elemTypes = mesh->elementTypes(); // global element types map<int, double> massFluxIntegral; // cellID -> integral double maxMassFluxIntegral = 0.0; double totalMassFlux = 0.0; double totalAbsMassFlux = 0.0; double maxCellMeasure = 0; double minCellMeasure = 1; for (vector< ElementTypePtr >::iterator elemTypeIt = elemTypes.begin(); elemTypeIt != elemTypes.end(); elemTypeIt++) { ElementTypePtr elemType = *elemTypeIt; vector< ElementPtr > elems = mesh->elementsOfTypeGlobal(elemType); vector<GlobalIndexType> cellIDs; for (int i=0; i<elems.size(); i++) { cellIDs.push_back(elems[i]->cellID()); } FieldContainer<double> physicalCellNodes = mesh->physicalCellNodesGlobal(elemType); BasisCachePtr basisCache = Teuchos::rcp( new BasisCache(elemType,mesh,polyOrder) ); // enrich by trial space order basisCache->setPhysicalCellNodes(physicalCellNodes,cellIDs,true); // true: create side caches FieldContainer<double> cellMeasures = basisCache->getCellMeasures(); FieldContainer<double> fakeRHSIntegrals(elems.size(),testOrdering->totalDofs()); massFluxTerm->integrate(fakeRHSIntegrals,testOrdering,basisCache,true); // true: force side evaluation // cout << "fakeRHSIntegrals:\n" << fakeRHSIntegrals; for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; // pick out the ones for testOne: massFluxIntegral[cellID] = fakeRHSIntegrals(i,testOneIndex); } // find the largest: for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; maxMassFluxIntegral = max(abs(massFluxIntegral[cellID]), maxMassFluxIntegral); } for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; maxCellMeasure = max(maxCellMeasure,cellMeasures(i)); minCellMeasure = min(minCellMeasure,cellMeasures(i)); maxMassFluxIntegral = max(abs(massFluxIntegral[cellID]), maxMassFluxIntegral); totalMassFlux += massFluxIntegral[cellID]; totalAbsMassFlux += abs( massFluxIntegral[cellID] ); } } if (rank==0) { cout << "largest mass flux: " << maxMassFluxIntegral << endl; cout << "total mass flux: " << totalMassFlux << endl; cout << "sum of mass flux absolute value: " << totalAbsMassFlux << endl; cout << "largest h: " << sqrt(maxCellMeasure) << endl; cout << "smallest h: " << sqrt(minCellMeasure) << endl; cout << "ratio of largest / smallest h: " << sqrt(maxCellMeasure) / sqrt(minCellMeasure) << endl; } if (rank == 0) { cout << "phi ID: " << phi->ID() << endl; cout << "psi1 ID: " << psi_1->ID() << endl; cout << "psi2 ID: " << psi_2->ID() << endl; cout << "streamMesh has " << streamMesh->numActiveElements() << " elements.\n"; cout << "solving for approximate stream function...\n"; } streamSolution->solve(useMumps); energyErrorTotal = streamSolution->energyErrorTotal(); if (rank == 0) { cout << "...solved.\n"; cout << "Stream mesh has energy error: " << energyErrorTotal << endl; } if (rank==0) { solution->writeToVTK("nsCavitySoln.vtk"); if (! meshHasTriangles ) { massFlux->writeBoundaryValuesToMATLABFile(solution->mesh(), "massFlux.dat"); u_mag->writeValuesToMATLABFile(solution->mesh(), "u_mag.m"); u_div->writeValuesToMATLABFile(solution->mesh(), "u_div.m"); solution->writeFieldsToFile(u1->ID(), "u1.m"); solution->writeFluxesToFile(u1hat->ID(), "u1_hat.dat"); solution->writeFieldsToFile(u2->ID(), "u2.m"); solution->writeFluxesToFile(u2hat->ID(), "u2_hat.dat"); solution->writeFieldsToFile(p->ID(), "p.m"); streamSolution->writeFieldsToFile(phi->ID(), "phi.m"); streamSolution->writeFluxesToFile(phi_hat->ID(), "phi_hat.dat"); streamSolution->writeFieldsToFile(psi_1->ID(), "psi1.m"); streamSolution->writeFieldsToFile(psi_2->ID(), "psi2.m"); vorticity->writeValuesToMATLABFile(streamMesh, "vorticity.m"); FunctionPtr ten = Teuchos::rcp( new ConstantScalarFunction(10) ); ten->writeBoundaryValuesToMATLABFile(solution->mesh(), "skeleton.dat"); cout << "wrote files: u_mag.m, u_div.m, u1.m, u1_hat.dat, u2.m, u2_hat.dat, p.m, phi.m, vorticity.m.\n"; } else { solution->writeToFile(u1->ID(), "u1.dat"); solution->writeToFile(u2->ID(), "u2.dat"); solution->writeToFile(u2->ID(), "p.dat"); cout << "wrote files: u1.dat, u2.dat, p.dat\n"; } FieldContainer<double> points = pointGrid(0, 1, 0, 1, 100); FieldContainer<double> pointData = solutionData(points, streamSolution, phi); GnuPlotUtil::writeXYPoints("phi_patch_navierStokes_cavity.dat", pointData); set<double> patchContourLevels = diagonalContourLevels(pointData,1); vector<string> patchDataPath; patchDataPath.push_back("phi_patch_navierStokes_cavity.dat"); GnuPlotUtil::writeContourPlotScript(patchContourLevels, patchDataPath, "lidCavityNavierStokes.p"); GnuPlotUtil::writeExactMeshSkeleton("lid_navierStokes_continuation_adaptive", mesh, 2); writePatchValues(0, 1, 0, 1, streamSolution, phi, "phi_patch.m"); writePatchValues(0, .1, 0, .1, streamSolution, phi, "phi_patch_detail.m"); writePatchValues(0, .01, 0, .01, streamSolution, phi, "phi_patch_minute_detail.m"); writePatchValues(0, .001, 0, .001, streamSolution, phi, "phi_patch_minute_minute_detail.m"); } return 0; }
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; }
void ExactSolution::L2NormOfError(FieldContainer<double> &errorSquaredPerCell, Solution &solution, ElementTypePtr elemTypePtr, int trialID, int sideIndex, int cubDegree, double solutionLift) { // BasisCache(ElementTypePtr elemType, Teuchos::RCP<Mesh> mesh = Teuchos::rcp( (Mesh*) NULL ), bool testVsTest=false, int cubatureDegreeEnrichment = 0) DofOrdering dofOrdering = *(elemTypePtr->trialOrderPtr.get()); BasisPtr basis = dofOrdering.getBasis(trialID,sideIndex); bool boundaryIntegral = solution.mesh()->bilinearForm()->isFluxOrTrace(trialID); BasisCachePtr basisCache; if (cubDegree <= 0) { // then take the default cub. degree basisCache = Teuchos::rcp( new BasisCache( elemTypePtr, solution.mesh() ) ); } else { // we could eliminate the logic below if we just added BasisCache::setCubatureDegree()... // (the logic below is just to make the enriched cubature match the requested cubature degree...) int maxTrialDegree; if (!boundaryIntegral) { maxTrialDegree = elemTypePtr->trialOrderPtr->maxBasisDegreeForVolume(); } else { maxTrialDegree = elemTypePtr->trialOrderPtr->maxBasisDegree(); // generally, this will be the trace degree } int maxTestDegree = elemTypePtr->testOrderPtr->maxBasisDegree(); int cubDegreeEnrichment = max(cubDegree - (maxTrialDegree + maxTestDegree), 0); basisCache = Teuchos::rcp( new BasisCache( elemTypePtr, solution.mesh(), false, cubDegreeEnrichment) ); } // much of this code is the same as what's in the volume integration in computeStiffness... FieldContainer<double> physicalCellNodes = solution.mesh()->physicalCellNodes(elemTypePtr); vector<GlobalIndexType> cellIDs = solution.mesh()->cellIDsOfType(elemTypePtr); basisCache->setPhysicalCellNodes(physicalCellNodes, cellIDs, true); if (boundaryIntegral) { basisCache = basisCache->getSideBasisCache(sideIndex); } FieldContainer<double> weightedMeasure = basisCache->getWeightedMeasures(); FieldContainer<double> weightedErrorSquared; int numCells = basisCache->getPhysicalCubaturePoints().dimension(0); int numCubPoints = basisCache->getPhysicalCubaturePoints().dimension(1); int spaceDim = basisCache->getPhysicalCubaturePoints().dimension(2); Teuchos::Array<int> dimensions; dimensions.push_back(numCells); dimensions.push_back(numCubPoints); int basisRank = BasisFactory::basisFactory()->getBasisRank(basis); if (basisRank==1) { dimensions.push_back(spaceDim); } FieldContainer<double> computedValues(dimensions); FieldContainer<double> exactValues(dimensions); if (solutionLift != 0.0) { int size = computedValues.size(); for (int i=0; i<size; i++) { computedValues[i] += solutionLift; } } solution.solutionValues(computedValues, trialID, basisCache); this->solutionValues(exactValues, trialID, basisCache); // cout << "ExactSolution: exact values:\n" << exactValues; // cout << "ExactSolution: computed values:\n" << computedValues; FieldContainer<double> errorSquared(numCells,numCubPoints); squaredDifference(errorSquared,computedValues,exactValues); weightedErrorSquared.resize(numCells,numCubPoints); for (int cellIndex=0; cellIndex<numCells; cellIndex++) { for (int ptIndex=0; ptIndex<numCubPoints; ptIndex++) { // following two lines for viewing in the debugger: double weight = weightedMeasure(cellIndex,ptIndex); double errorSquaredVal = errorSquared(cellIndex,ptIndex); weightedErrorSquared(cellIndex,ptIndex) = errorSquared(cellIndex,ptIndex) * weightedMeasure(cellIndex,ptIndex); } } // compute the integral errorSquaredPerCell.initialize(0.0); int numPoints = weightedErrorSquared.dimension(1); for (int cellIndex=0; cellIndex<numCells; cellIndex++) { for (int ptIndex=0; ptIndex<numPoints; ptIndex++) { errorSquaredPerCell(cellIndex) += weightedErrorSquared(cellIndex,ptIndex); } } }
bool LobattoBasisTests::testSimpleStiffnessMatrix() { bool success = true; int rank = Teuchos::GlobalMPISession::getRank(); VarFactory varFactory; VarPtr u = varFactory.fieldVar("u"); VarPtr un = varFactory.fluxVar("un_hat"); VarPtr v = varFactory.testVar("v", HGRAD); BFPtr bf = Teuchos::rcp( new BF(varFactory) ); vector<double> beta; beta.push_back(1); beta.push_back(1); bf->addTerm(beta * u, v->grad()); bf->addTerm(un, v); DofOrderingPtr trialSpace = Teuchos::rcp( new DofOrdering ); DofOrderingPtr testSpace = Teuchos::rcp( new DofOrdering ); const int numSides = 4; const int spaceDim = 2; int fieldOrder = 3; int testOrder = fieldOrder+2; BasisPtr fieldBasis = Camellia::intrepidQuadHGRAD(fieldOrder); BasisPtr fluxBasis = Camellia::intrepidLineHGRAD(fieldOrder); trialSpace->addEntry(u->ID(), fieldBasis, fieldBasis->rangeRank()); for (int i=0; i<numSides; i++) { trialSpace->addEntry(un->ID(), fluxBasis, fluxBasis->rangeRank(), i); } BasisPtr testBasis = Camellia::lobattoQuadHGRAD(testOrder+1,false); // +1 because it lives in HGRAD testSpace->addEntry(v->ID(), testBasis, testBasis->rangeRank()); int numTrialDofs = trialSpace->totalDofs(); int numTestDofs = testSpace->totalDofs(); int numCells = 1; FieldContainer<double> cellNodes(numCells,numSides,spaceDim); cellNodes(0,0,0) = 0; cellNodes(0,0,1) = 0; cellNodes(0,1,0) = 1; cellNodes(0,1,1) = 0; cellNodes(0,2,0) = 1; cellNodes(0,2,1) = 1; cellNodes(0,3,0) = 0; cellNodes(0,3,1) = 1; FieldContainer<double> stiffness(numCells,numTestDofs,numTrialDofs); FieldContainer<double> cellSideParities(numCells,numSides); cellSideParities.initialize(1.0); Teuchos::RCP<shards::CellTopology> quad_4 = Teuchos::rcp( new shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() ) ); Teuchos::RCP<ElementType> elemType = Teuchos::rcp( new ElementType(trialSpace, testSpace, quad_4)); BasisCachePtr basisCache = Teuchos::rcp( new BasisCache(elemType) ); vector<GlobalIndexType> cellIDs; cellIDs.push_back(0); basisCache->setPhysicalCellNodes(cellNodes, cellIDs, true); bf->stiffnessMatrix(stiffness, elemType, cellSideParities, basisCache); // TODO: finish this test // cout << stiffness; if (rank==0) cout << "Warning: testSimpleStiffnessMatrix() unfinished.\n"; return success; }
bool MeshTestUtility::determineRefTestPointsForNeighbors(MeshTopologyViewPtr meshTopo, CellPtr fineCell, unsigned int sideOrdinal, FieldContainer<double> &fineSideRefPoints, FieldContainer<double> &fineCellRefPoints, FieldContainer<double> &coarseSideRefPoints, FieldContainer<double> &coarseCellRefPoints) { unsigned spaceDim = meshTopo->getDimension(); unsigned sideDim = spaceDim - 1; if (spaceDim == 1) { fineSideRefPoints.resize(0,0); coarseSideRefPoints.resize(0,0); fineCellRefPoints.resize(1,1); coarseCellRefPoints.resize(1,1); FieldContainer<double> lineRefNodes(2,1); CellTopoPtr line = CellTopology::line(); CamelliaCellTools::refCellNodesForTopology(lineRefNodes, line); fineCellRefPoints[0] = lineRefNodes[sideOrdinal]; unsigned neighborSideOrdinal = fineCell->getNeighborInfo(sideOrdinal, meshTopo).second; if (neighborSideOrdinal != -1) { coarseCellRefPoints[0] = lineRefNodes[neighborSideOrdinal]; return true; } else { return false; } } pair<GlobalIndexType, unsigned> neighborInfo = fineCell->getNeighborInfo(sideOrdinal, meshTopo); if (neighborInfo.first == -1) { // boundary return false; } CellPtr neighborCell = meshTopo->getCell(neighborInfo.first); if (neighborCell->isParent(meshTopo)) { return false; // fineCell isn't the finer of the two... } CellTopoPtr fineSideTopo = fineCell->topology()->getSubcell(sideDim, sideOrdinal); CubatureFactory cubFactory; int cubDegree = 4; // fairly arbitrary choice: enough to get a decent number of test points... Teuchos::RCP<Cubature<double> > fineSideTopoCub = cubFactory.create(fineSideTopo, cubDegree); int numCubPoints = fineSideTopoCub->getNumPoints(); FieldContainer<double> cubPoints(numCubPoints, sideDim); FieldContainer<double> cubWeights(numCubPoints); // we neglect these... fineSideTopoCub->getCubature(cubPoints, cubWeights); FieldContainer<double> sideRefCellNodes(fineSideTopo->getNodeCount(),sideDim); CamelliaCellTools::refCellNodesForTopology(sideRefCellNodes, fineSideTopo); int numTestPoints = numCubPoints + fineSideTopo->getNodeCount(); FieldContainer<double> testPoints(numTestPoints, sideDim); for (int ptOrdinal=0; ptOrdinal<testPoints.dimension(0); ptOrdinal++) { if (ptOrdinal<fineSideTopo->getNodeCount()) { for (int d=0; d<sideDim; d++) { testPoints(ptOrdinal,d) = sideRefCellNodes(ptOrdinal,d); } } else { for (int d=0; d<sideDim; d++) { testPoints(ptOrdinal,d) = cubPoints(ptOrdinal-fineSideTopo->getNodeCount(),d); } } } fineSideRefPoints = testPoints; fineCellRefPoints.resize(numTestPoints, spaceDim); CamelliaCellTools::mapToReferenceSubcell(fineCellRefPoints, testPoints, sideDim, sideOrdinal, fineCell->topology()); CellTopoPtr coarseSideTopo = neighborCell->topology()->getSubcell(sideDim, neighborInfo.second); unsigned fineSideAncestorPermutation = fineCell->ancestralPermutationForSubcell(sideDim, sideOrdinal, meshTopo); unsigned coarseSidePermutation = neighborCell->subcellPermutation(sideDim, neighborInfo.second); unsigned coarseSideAncestorPermutationInverse = CamelliaCellTools::permutationInverse(coarseSideTopo, coarseSidePermutation); unsigned composedPermutation = CamelliaCellTools::permutationComposition(coarseSideTopo, coarseSideAncestorPermutationInverse, fineSideAncestorPermutation); // goes from coarse ordering to fine. RefinementBranch fineRefBranch = fineCell->refinementBranchForSide(sideOrdinal, meshTopo); FieldContainer<double> fineSideNodes(fineSideTopo->getNodeCount(), sideDim); // relative to the ancestral cell whose neighbor is compatible if (fineRefBranch.size() == 0) { CamelliaCellTools::refCellNodesForTopology(fineSideNodes, coarseSideTopo, composedPermutation); } else { FieldContainer<double> ancestralSideNodes(coarseSideTopo->getNodeCount(), sideDim); CamelliaCellTools::refCellNodesForTopology(ancestralSideNodes, coarseSideTopo, composedPermutation); RefinementBranch fineSideRefBranch = RefinementPattern::sideRefinementBranch(fineRefBranch, sideOrdinal); fineSideNodes = RefinementPattern::descendantNodes(fineSideRefBranch, ancestralSideNodes); } BasisCachePtr sideTopoCache = Teuchos::rcp( new BasisCache(fineSideTopo, 1, false) ); sideTopoCache->setRefCellPoints(testPoints); // add cell dimension fineSideNodes.resize(1,fineSideNodes.dimension(0), fineSideNodes.dimension(1)); sideTopoCache->setPhysicalCellNodes(fineSideNodes, vector<GlobalIndexType>(), false); coarseSideRefPoints = sideTopoCache->getPhysicalCubaturePoints(); // strip off the cell dimension: coarseSideRefPoints.resize(coarseSideRefPoints.dimension(1),coarseSideRefPoints.dimension(2)); coarseCellRefPoints.resize(numTestPoints,spaceDim); CamelliaCellTools::mapToReferenceSubcell(coarseCellRefPoints, coarseSideRefPoints, sideDim, neighborInfo.second, neighborCell->topology()); return true; // containers filled.... }
bool HConvergenceStudyTests::testBestApproximationErrorComputation() { bool success = true; bool enrichVelocity = false; // true would be for the "compliant" norm, which isn't working well yet int minLogElements = 0, maxLogElements = minLogElements; int numCells1D = pow(2.0,minLogElements); int H1Order = 1; int pToAdd = 2; double tol = 1e-16; double Re = 40.0; VarFactory varFactory = VGPStokesFormulation::vgpVarFactory(); VarPtr u1_vgp = varFactory.fieldVar(VGP_U1_S); VarPtr u2_vgp = varFactory.fieldVar(VGP_U2_S); VarPtr sigma11_vgp = varFactory.fieldVar(VGP_SIGMA11_S); VarPtr sigma12_vgp = varFactory.fieldVar(VGP_SIGMA12_S); VarPtr sigma21_vgp = varFactory.fieldVar(VGP_SIGMA21_S); VarPtr sigma22_vgp = varFactory.fieldVar(VGP_SIGMA22_S); VarPtr p_vgp = varFactory.fieldVar(VGP_P_S); VGPStokesFormulation stokesForm(1/Re); int numCellsFineMesh = 20; // for computing a zero-mean pressure int H1OrderFineMesh = 5; // define Kovasznay domain: FieldContainer<double> quadPointsKovasznay(4,2); // Domain from Evans Hughes for Navier-Stokes: quadPointsKovasznay(0,0) = 0.0; // x1 quadPointsKovasznay(0,1) = -0.5; // y1 quadPointsKovasznay(1,0) = 1.0; quadPointsKovasznay(1,1) = -0.5; quadPointsKovasznay(2,0) = 1.0; quadPointsKovasznay(2,1) = 0.5; quadPointsKovasznay(3,0) = 0.0; quadPointsKovasznay(3,1) = 0.5; FunctionPtr zero = Function::zero(); bool dontEnhanceFluxes = false; VGPNavierStokesProblem zeroProblem = VGPNavierStokesProblem(Re, quadPointsKovasznay, numCellsFineMesh, numCellsFineMesh, H1OrderFineMesh, pToAdd, zero, zero, zero, enrichVelocity, dontEnhanceFluxes); FunctionPtr u1_exact, u2_exact, p_exact; NavierStokesFormulation::setKovasznay(Re, zeroProblem.mesh(), u1_exact, u2_exact, p_exact); VGPNavierStokesProblem problem = VGPNavierStokesProblem(Re,quadPointsKovasznay, numCells1D,numCells1D, H1Order, pToAdd, u1_exact, u2_exact, p_exact, enrichVelocity, dontEnhanceFluxes); HConvergenceStudy study(problem.exactSolution(), problem.mesh()->bilinearForm(), problem.exactSolution()->rhs(), problem.backgroundFlow()->bc(), problem.bf()->graphNorm(), minLogElements, maxLogElements, H1Order, pToAdd, false, false, false); study.setReportRelativeErrors(false); // we want absolute errors Teuchos::RCP<Mesh> mesh = problem.mesh(); int cubatureDegreeEnrichment = 10; int L2Order = H1Order - 1; int meshCubatureDegree = L2Order + H1Order + pToAdd; study.setCubatureDegreeForExact(cubatureDegreeEnrichment + meshCubatureDegree); FunctionPtr f = u1_exact; int trialID = u1_vgp->ID(); { double fIntegral = f->integrate(mesh,cubatureDegreeEnrichment); // cout << "testBestApproximationErrorComputation: integral of f on whole mesh = " << fIntegral << endl; double l2ErrorOfAverage = (Function::constant(fIntegral) - f)->l2norm(mesh,cubatureDegreeEnrichment); // cout << "testBestApproximationErrorComputation: l2 error of fIntegral: " << l2ErrorOfAverage << endl; ElementTypePtr elemType = mesh->elementTypes()[0]; vector<GlobalIndexType> cellIDs = mesh->cellIDsOfTypeGlobal(elemType); bool testVsTest = false; BasisCachePtr basisCache = Teuchos::rcp( new BasisCache(elemType, mesh, testVsTest, cubatureDegreeEnrichment) ); basisCache->setPhysicalCellNodes(mesh->physicalCellNodesGlobal(elemType), cellIDs, false); // false: no side cache FieldContainer<double> projectionValues(cellIDs.size()); f->integrate(projectionValues, basisCache); FieldContainer<double> cellMeasures = basisCache->getCellMeasures(); for (int i=0; i<projectionValues.size(); i++) { projectionValues(i) /= cellMeasures(i); } // since we're not worried about the actual solution values at all, just use a single zero solution: vector< SolutionPtr > solutions; solutions.push_back( problem.backgroundFlow() ); study.setSolutions(solutions); // this will call computeError() double approximationError = study.bestApproximationErrors()[trialID][0]; // 0: solution/mesh index // for a single-cell mesh, approximation error should be the same as the L^2 error of the average double diff = abs(approximationError - l2ErrorOfAverage); if (diff > tol) { cout << "testBestApproximationErrorComputation: diff " << diff << " exceeds tol " << tol << endl; success = false; } else { // cout << "testBestApproximationErrorComputation: diff " << diff << " is below tol " << tol << endl; } } return success; }
bool LinearTermTests::testRieszInversionAsProjection() { bool success = true; //////////////////// DECLARE VARIABLES /////////////////////// // define test variables VarFactoryPtr varFactory = VarFactory::varFactory(); VarPtr tau = varFactory->testVar("\\tau", HDIV); VarPtr v = varFactory->testVar("v", HGRAD); // define trial variables VarPtr uhat = varFactory->traceVar("\\widehat{u}"); VarPtr beta_n_u_minus_sigma_n = varFactory->fluxVar("\\widehat{\\beta \\cdot n u - \\sigma_{n}}"); VarPtr u = varFactory->fieldVar("u"); VarPtr sigma1 = varFactory->fieldVar("\\sigma_1"); VarPtr sigma2 = varFactory->fieldVar("\\sigma_2"); vector<double> beta; beta.push_back(1.0); beta.push_back(0.0); double eps = .01; //////////////////// DEFINE BILINEAR FORM /////////////////////// BFPtr confusionBF = Teuchos::rcp( new BF(varFactory) ); // tau terms: confusionBF->addTerm(sigma1 / eps, tau->x()); confusionBF->addTerm(sigma2 / eps, tau->y()); confusionBF->addTerm(u, tau->div()); confusionBF->addTerm(uhat, -tau->dot_normal()); // v terms: confusionBF->addTerm( sigma1, v->dx() ); confusionBF->addTerm( sigma2, v->dy() ); confusionBF->addTerm( -u, beta * v->grad() ); confusionBF->addTerm( beta_n_u_minus_sigma_n, v); //////////////////// BUILD MESH /////////////////////// // define nodes for mesh int H1Order = 2; int pToAdd = 2; 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 nCells = 2; int horizontalCells = nCells, verticalCells = nCells; // create a new mesh: MeshPtr myMesh = MeshFactory::buildQuadMesh(quadPoints, horizontalCells, verticalCells, confusionBF, H1Order, H1Order+pToAdd); ElementTypePtr elemType = myMesh->getElement(0)->elementType(); BasisCachePtr basisCache = Teuchos::rcp(new BasisCache(elemType, myMesh)); vector<GlobalIndexType> cellIDs = myMesh->cellIDsOfTypeGlobal(elemType); bool createSideCacheToo = true; basisCache->setPhysicalCellNodes(myMesh->physicalCellNodesGlobal(elemType), cellIDs, createSideCacheToo); LinearTermPtr integrand = Teuchos::rcp(new LinearTerm); // residual FunctionPtr x = Function::xn(1); FunctionPtr y = Function::yn(1); FunctionPtr testFxn1 = x; FunctionPtr testFxn2 = y; FunctionPtr fxnToProject = x * y + 1.0; integrand->addTerm(fxnToProject * v); IPPtr ip_L2 = Teuchos::rcp(new IP); ip_L2->addTerm(v); ip_L2->addTerm(tau); Teuchos::RCP<RieszRep> riesz = Teuchos::rcp(new RieszRep(myMesh, ip_L2, integrand)); riesz->computeRieszRep(); FunctionPtr rieszFxn = RieszRep::repFunction(v,riesz); int numCells = basisCache->getPhysicalCubaturePoints().dimension(0); int numPts = basisCache->getPhysicalCubaturePoints().dimension(1); FieldContainer<double> valProject( numCells, numPts ); FieldContainer<double> valExpected( numCells, numPts ); rieszFxn->values(valProject,basisCache); fxnToProject->values(valExpected,basisCache); // int rank = Teuchos::GlobalMPISession::getRank(); // if (rank==0) cout << "physicalCubaturePoints:\n" << basisCache->getPhysicalCubaturePoints(); double maxDiff; double tol = 1e-12; success = TestSuite::fcsAgree(valProject,valExpected,tol,maxDiff); if (success==false) { cout << "Failed Riesz Inversion Projection test with maxDiff = " << maxDiff << endl; serializeOutput("valExpected", valExpected); serializeOutput("valProject", valProject); serializeOutput("physicalPoints", basisCache->getPhysicalCubaturePoints()); } return allSuccess(success); }
int main(int argc, char *argv[]) { // Process command line arguments if (argc > 1) numRefs = atof(argv[1]); #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 FunctionPtr beta = Teuchos::rcp(new Beta()); //////////////////////////////////////////////////////////////////// // DEFINE VARIABLES //////////////////////////////////////////////////////////////////// // test variables VarFactory varFactory; VarPtr tau = varFactory.testVar("\\tau", HDIV); VarPtr v = varFactory.testVar("v", HGRAD); // trial variables VarPtr uhat = varFactory.traceVar("\\widehat{u}"); VarPtr beta_n_u_minus_sigma_n = varFactory.fluxVar("\\widehat{\\beta \\cdot n u - \\sigma_{n}}"); VarPtr u = varFactory.fieldVar("u"); VarPtr sigma1 = varFactory.fieldVar("\\sigma_1"); VarPtr sigma2 = varFactory.fieldVar("\\sigma_2"); //////////////////////////////////////////////////////////////////// // CREATE MESH //////////////////////////////////////////////////////////////////// BFPtr confusionBF = Teuchos::rcp( new BF(varFactory) ); FieldContainer<double> meshBoundary(4,2); meshBoundary(0,0) = 0.0; // x1 meshBoundary(0,1) = -2.0; // y1 meshBoundary(1,0) = 4.0; meshBoundary(1,1) = -2.0; meshBoundary(2,0) = 4.0; meshBoundary(2,1) = 2.0; meshBoundary(3,0) = 0.0; meshBoundary(3,1) = 2.0; int horizontalCells = 4, verticalCells = 4; // create a pointer to a new mesh: Teuchos::RCP<Mesh> mesh = Mesh::buildQuadMesh(meshBoundary, horizontalCells, verticalCells, confusionBF, H1Order, H1Order+pToAdd, false); //////////////////////////////////////////////////////////////////// // INITIALIZE BACKGROUND FLOW FUNCTIONS //////////////////////////////////////////////////////////////////// BCPtr nullBC = Teuchos::rcp((BC*)NULL); RHSPtr nullRHS = Teuchos::rcp((RHS*)NULL); IPPtr nullIP = Teuchos::rcp((IP*)NULL); SolutionPtr prevTimeFlow = Teuchos::rcp(new Solution(mesh, nullBC, nullRHS, nullIP) ); SolutionPtr flowResidual = Teuchos::rcp(new Solution(mesh, nullBC, nullRHS, nullIP) ); FunctionPtr u_prev_time = Teuchos::rcp( new PreviousSolutionFunction(prevTimeFlow, u) ); // ==================== SET INITIAL GUESS ========================== double u_free = 0.0; double sigma1_free = 0.0; double sigma2_free = 0.0; map<int, Teuchos::RCP<Function> > functionMap; functionMap[u->ID()] = Teuchos::rcp( new ConstantScalarFunction(u_free) ); functionMap[sigma1->ID()] = Teuchos::rcp( new ConstantScalarFunction(sigma1_free) ); functionMap[sigma2->ID()] = Teuchos::rcp( new ConstantScalarFunction(sigma2_free) ); prevTimeFlow->projectOntoMesh(functionMap); // ==================== END SET INITIAL GUESS ========================== //////////////////////////////////////////////////////////////////// // DEFINE BILINEAR FORM //////////////////////////////////////////////////////////////////// // tau terms: confusionBF->addTerm(sigma1 / epsilon, tau->x()); confusionBF->addTerm(sigma2 / epsilon, tau->y()); confusionBF->addTerm(u, tau->div()); confusionBF->addTerm(-uhat, tau->dot_normal()); // v terms: confusionBF->addTerm( sigma1, v->dx() ); confusionBF->addTerm( sigma2, v->dy() ); confusionBF->addTerm( beta * u, - v->grad() ); confusionBF->addTerm( beta_n_u_minus_sigma_n, v); //////////////////////////////////////////////////////////////////// // TIMESTEPPING TERMS //////////////////////////////////////////////////////////////////// Teuchos::RCP<RHSEasy> rhs = Teuchos::rcp( new RHSEasy ); double dt = 0.25; FunctionPtr invDt = Teuchos::rcp(new ScalarParamFunction(1.0/dt)); if (rank==0){ cout << "Timestep dt = " << dt << endl; } if (transient) { confusionBF->addTerm( u, invDt*v ); rhs->addTerm( u_prev_time * invDt * v ); } //////////////////////////////////////////////////////////////////// // DEFINE INNER PRODUCT //////////////////////////////////////////////////////////////////// // mathematician's norm IPPtr mathIP = Teuchos::rcp(new IP()); mathIP->addTerm(tau); mathIP->addTerm(tau->div()); mathIP->addTerm(v); mathIP->addTerm(v->grad()); // quasi-optimal norm IPPtr qoptIP = Teuchos::rcp(new IP); qoptIP->addTerm( v ); qoptIP->addTerm( tau / epsilon + v->grad() ); qoptIP->addTerm( beta * v->grad() - tau->div() ); // robust test norm IPPtr robIP = Teuchos::rcp(new IP); FunctionPtr ip_scaling = Teuchos::rcp( new EpsilonScaling(epsilon) ); if (!enforceLocalConservation) { robIP->addTerm( ip_scaling * v ); if (transient) robIP->addTerm( invDt * v ); } robIP->addTerm( sqrt(epsilon) * v->grad() ); // Weight these two terms for inflow FunctionPtr ip_weight = Teuchos::rcp( new IPWeight() ); robIP->addTerm( ip_weight * beta * v->grad() ); robIP->addTerm( ip_weight * tau->div() ); robIP->addTerm( ip_scaling/sqrt(epsilon) * tau ); if (enforceLocalConservation) robIP->addZeroMeanTerm( v ); //////////////////////////////////////////////////////////////////// // DEFINE RHS //////////////////////////////////////////////////////////////////// FunctionPtr f = Teuchos::rcp( new ConstantScalarFunction(0.0) ); rhs->addTerm( f * v ); // obviously, with f = 0 adding this term is not necessary! //////////////////////////////////////////////////////////////////// // DEFINE BC //////////////////////////////////////////////////////////////////// Teuchos::RCP<BCEasy> bc = Teuchos::rcp( new BCEasy ); // Teuchos::RCP<PenaltyConstraints> pc = Teuchos::rcp( new PenaltyConstraints ); SpatialFilterPtr lBoundary = Teuchos::rcp( new LeftBoundary ); SpatialFilterPtr tbBoundary = Teuchos::rcp( new TopBottomBoundary ); SpatialFilterPtr rBoundary = Teuchos::rcp( new RightBoundary ); FunctionPtr u0 = Teuchos::rcp( new ZeroBC ); FunctionPtr u_inlet = Teuchos::rcp( new InletBC ); // FunctionPtr n = Teuchos::rcp( new UnitNormalFunction ); bc->addDirichlet(beta_n_u_minus_sigma_n, lBoundary, u_inlet); bc->addDirichlet(beta_n_u_minus_sigma_n, tbBoundary, u0); bc->addDirichlet(uhat, rBoundary, u0); // pc->addConstraint(beta_n_u_minus_sigma_n - uhat == u0, rBoundary); //////////////////////////////////////////////////////////////////// // CREATE SOLUTION OBJECT //////////////////////////////////////////////////////////////////// Teuchos::RCP<Solution> solution = Teuchos::rcp( new Solution(mesh, bc, rhs, robIP) ); // solution->setFilter(pc); // ==================== Enforce Local Conservation ================== if (enforceLocalConservation) { if (transient) { FunctionPtr conserved_rhs = u_prev_time * invDt; LinearTermPtr conserved_quantity = invDt * u; LinearTermPtr flux_part = Teuchos::rcp(new LinearTerm(-1.0, beta_n_u_minus_sigma_n)); conserved_quantity->addTerm(flux_part, true); // conserved_quantity = conserved_quantity - beta_n_u_minus_sigma_n; solution->lagrangeConstraints()->addConstraint(conserved_quantity == conserved_rhs); } else { FunctionPtr zero = Teuchos::rcp( new ConstantScalarFunction(0.0) ); solution->lagrangeConstraints()->addConstraint(beta_n_u_minus_sigma_n == zero); } } // ==================== Register Solutions ========================== mesh->registerSolution(solution); mesh->registerSolution(prevTimeFlow); // u_t(i-1) mesh->registerSolution(flowResidual); // u_t(i-1) double energyThreshold = 0.25; // for mesh refinements Teuchos::RCP<RefinementStrategy> refinementStrategy; refinementStrategy = Teuchos::rcp(new RefinementStrategy(solution,energyThreshold)); //////////////////////////////////////////////////////////////////// // PSEUDO-TIME SOLVE STRATEGY //////////////////////////////////////////////////////////////////// double time_tol = 1e-8; for (int refIndex=0; refIndex<=numRefs; refIndex++) { double L2_time_residual = 1e7; int timestepCount = 0; if (!transient) numTimeSteps = 1; while((L2_time_residual > time_tol) && (timestepCount < numTimeSteps)) { solution->solve(false); // subtract solutions to get residual flowResidual->setSolution(solution); // reset previous time solution to current time sol flowResidual->addSolution(prevTimeFlow, -1.0); double L2u = flowResidual->L2NormOfSolutionGlobal(u->ID()); double L2sigma1 = flowResidual->L2NormOfSolutionGlobal(sigma1->ID()); double L2sigma2 = flowResidual->L2NormOfSolutionGlobal(sigma2->ID()); L2_time_residual = sqrt(L2u*L2u + L2sigma1*L2sigma1 + L2sigma2*L2sigma2); cout << endl << "Timestep: " << timestepCount << ", dt = " << dt << ", Time residual = " << L2_time_residual << endl; if (rank == 0) { stringstream outfile; if (transient) outfile << "TransientConfusion_" << refIndex << "_" << timestepCount; else outfile << "TransientConfusion_" << refIndex; solution->writeToVTK(outfile.str(), 5); } ////////////////////////////////////////////////////////////////////////// // Check conservation by testing against one ////////////////////////////////////////////////////////////////////////// VarPtr testOne = varFactory.testVar("1", CONSTANT_SCALAR); // Create a fake bilinear form for the testing BFPtr fakeBF = Teuchos::rcp( new BF(varFactory) ); // Define our mass flux FunctionPtr flux_current_time = Teuchos::rcp( new PreviousSolutionFunction(solution, beta_n_u_minus_sigma_n) ); FunctionPtr delta_u = Teuchos::rcp( new PreviousSolutionFunction(flowResidual, u) ); LinearTermPtr surfaceFlux = -1.0 * flux_current_time * testOne; LinearTermPtr volumeChange = invDt * delta_u * testOne; LinearTermPtr massFluxTerm; if (transient) { massFluxTerm = volumeChange; // massFluxTerm->addTerm(surfaceFlux); } else { massFluxTerm = surfaceFlux; } // cout << "surface case = " << surfaceFlux->summands()[0].first->boundaryValueOnly() << " volume case = " << volumeChange->summands()[0].first->boundaryValueOnly() << endl; // FunctionPtr massFlux= Teuchos::rcp( new PreviousSolutionFunction(solution, beta_n_u_minus_sigma_n) ); // LinearTermPtr massFluxTerm = massFlux * testOne; Teuchos::RCP<shards::CellTopology> quadTopoPtr = Teuchos::rcp(new shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() )); DofOrderingFactory dofOrderingFactory(fakeBF); int fakeTestOrder = H1Order; DofOrderingPtr testOrdering = dofOrderingFactory.testOrdering(fakeTestOrder, *quadTopoPtr); int testOneIndex = testOrdering->getDofIndex(testOne->ID(),0); vector< ElementTypePtr > elemTypes = mesh->elementTypes(); // global element types map<int, double> massFluxIntegral; // cellID -> integral double maxMassFluxIntegral = 0.0; double totalMassFlux = 0.0; double totalAbsMassFlux = 0.0; for (vector< ElementTypePtr >::iterator elemTypeIt = elemTypes.begin(); elemTypeIt != elemTypes.end(); elemTypeIt++) { ElementTypePtr elemType = *elemTypeIt; vector< ElementPtr > elems = mesh->elementsOfTypeGlobal(elemType); vector<int> cellIDs; for (int i=0; i<elems.size(); i++) { cellIDs.push_back(elems[i]->cellID()); } FieldContainer<double> physicalCellNodes = mesh->physicalCellNodesGlobal(elemType); BasisCachePtr basisCache = Teuchos::rcp( new BasisCache(elemType,mesh) ); basisCache->setPhysicalCellNodes(physicalCellNodes,cellIDs,true); // true: create side caches FieldContainer<double> cellMeasures = basisCache->getCellMeasures(); FieldContainer<double> fakeRHSIntegrals(elems.size(),testOrdering->totalDofs()); massFluxTerm->integrate(fakeRHSIntegrals,testOrdering,basisCache,true); // true: force side evaluation for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; // pick out the ones for testOne: massFluxIntegral[cellID] = fakeRHSIntegrals(i,testOneIndex); } // find the largest: for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; maxMassFluxIntegral = max(abs(massFluxIntegral[cellID]), maxMassFluxIntegral); } for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; maxMassFluxIntegral = max(abs(massFluxIntegral[cellID]), maxMassFluxIntegral); totalMassFlux += massFluxIntegral[cellID]; totalAbsMassFlux += abs( massFluxIntegral[cellID] ); } } // Print results from processor with rank 0 if (rank == 0) { cout << "largest mass flux: " << maxMassFluxIntegral << endl; cout << "total mass flux: " << totalMassFlux << endl; cout << "sum of mass flux absolute value: " << totalAbsMassFlux << endl; } prevTimeFlow->setSolution(solution); // reset previous time solution to current time sol timestepCount++; } if (refIndex < numRefs){ if (rank==0){ cout << "Performing refinement number " << refIndex << endl; } refinementStrategy->refine(rank==0); // RESET solution every refinement - make sure discretization error doesn't creep in // prevTimeFlow->projectOntoMesh(functionMap); } } return 0; }
// tests Riesz inversion by integration by parts bool LinearTermTests::testRieszInversion() { bool success = true; //////////////////// DECLARE VARIABLES /////////////////////// // define test variables VarFactoryPtr varFactory = VarFactory::varFactory(); VarPtr tau = varFactory->testVar("\\tau", HDIV); VarPtr v = varFactory->testVar("v", HGRAD); // define trial variables VarPtr uhat = varFactory->traceVar("\\widehat{u}"); VarPtr beta_n_u_minus_sigma_n = varFactory->fluxVar("\\widehat{\\beta \\cdot n u - \\sigma_{n}}"); VarPtr u = varFactory->fieldVar("u"); VarPtr sigma1 = varFactory->fieldVar("\\sigma_1"); VarPtr sigma2 = varFactory->fieldVar("\\sigma_2"); vector<double> beta; beta.push_back(1.0); beta.push_back(0.0); double eps = .01; //////////////////// DEFINE BILINEAR FORM /////////////////////// BFPtr confusionBF = Teuchos::rcp( new BF(varFactory) ); // tau terms: confusionBF->addTerm(sigma1 / eps, tau->x()); confusionBF->addTerm(sigma2 / eps, tau->y()); confusionBF->addTerm(u, tau->div()); confusionBF->addTerm(uhat, -tau->dot_normal()); // v terms: confusionBF->addTerm( sigma1, v->dx() ); confusionBF->addTerm( sigma2, v->dy() ); confusionBF->addTerm( -u, beta * v->grad() ); confusionBF->addTerm( beta_n_u_minus_sigma_n, v); //////////////////// BUILD MESH /////////////////////// // define nodes for mesh int H1Order = 1; int pToAdd = 1; 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 nCells = 1; int horizontalCells = nCells, verticalCells = nCells; // create a pointer to a new mesh: Teuchos::RCP<Mesh> myMesh = MeshFactory::buildQuadMesh(quadPoints, horizontalCells, verticalCells, confusionBF, H1Order, H1Order+pToAdd); ElementTypePtr elemType = myMesh->getElement(0)->elementType(); BasisCachePtr basisCache = Teuchos::rcp(new BasisCache(elemType, myMesh)); vector<GlobalIndexType> cellIDs; vector<ElementPtr> elems = myMesh->activeElements(); vector<ElementPtr>::iterator elemIt; for (elemIt=elems.begin(); elemIt!=elems.end(); elemIt++) { int cellID = (*elemIt)->cellID(); cellIDs.push_back(cellID); } bool createSideCacheToo = true; basisCache->setPhysicalCellNodes(myMesh->physicalCellNodesGlobal(elemType), cellIDs, createSideCacheToo); LinearTermPtr integrand = Teuchos::rcp(new LinearTerm);// residual LinearTermPtr integrandIBP = Teuchos::rcp(new LinearTerm);// residual vector<double> e1(2); // (1,0) vector<double> e2(2); // (0,1) e1[0] = 1; e2[1] = 1; FunctionPtr n = Function::normal(); FunctionPtr X = Function::xn(1); FunctionPtr Y = Function::yn(1); FunctionPtr testFxn1 = X; FunctionPtr testFxn2 = Y; FunctionPtr divTestFxn = testFxn1->dx() + testFxn2->dy(); FunctionPtr vectorTest = testFxn1*e1 + testFxn2*e2; integrand->addTerm(divTestFxn*v); integrandIBP->addTerm(vectorTest*n*v - vectorTest*v->grad()); // boundary term IPPtr sobolevIP = Teuchos::rcp(new IP); sobolevIP->addTerm(v); sobolevIP->addTerm(tau); Teuchos::RCP<RieszRep> riesz = Teuchos::rcp(new RieszRep(myMesh, sobolevIP, integrand)); // riesz->setPrintOption(true); riesz->computeRieszRep(); Teuchos::RCP<RieszRep> rieszIBP = Teuchos::rcp(new RieszRep(myMesh, sobolevIP, integrandIBP)); riesz->setFunctional(integrandIBP); // rieszIBP->setPrintOption(true); rieszIBP->computeRieszRep(); FunctionPtr rieszOrigFxn = RieszRep::repFunction(v,riesz); FunctionPtr rieszIBPFxn = RieszRep::repFunction(v,rieszIBP); int numCells = basisCache->getPhysicalCubaturePoints().dimension(0); int numPts = basisCache->getPhysicalCubaturePoints().dimension(1); FieldContainer<double> valOriginal( numCells, numPts); FieldContainer<double> valIBP( numCells, numPts); rieszOrigFxn->values(valOriginal,basisCache); rieszIBPFxn->values(valIBP,basisCache); double maxDiff; double tol = 1e-14; success = TestSuite::fcsAgree(valOriginal,valIBP,tol,maxDiff); if (success==false) { cout << "Failed TestRieszInversion with maxDiff = " << maxDiff << endl; } return success; }
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: bool useTriangles = false; FieldContainer<double> meshPoints(4,2); meshPoints(0,0) = 0.0; // x1 meshPoints(0,1) = 0.0; // y1 meshPoints(1,0) = 1.0; meshPoints(1,1) = 0.0; meshPoints(2,0) = 1.0; meshPoints(2,1) = 1.0; meshPoints(3,0) = 0.0; meshPoints(3,1) = 1.0; int H1Order = polyOrder + 1; int horizontalCells = 4, verticalCells = 4; double energyThreshold = 0.2; // for mesh refinements double nonlinearStepSize = 0.5; double nonlinearRelativeEnergyTolerance = 1e-8; // used to determine convergence of the nonlinear solution //////////////////////////////////////////////////////////////////// // DEFINE VARIABLES //////////////////////////////////////////////////////////////////// // new-style bilinear form definition VarFactory varFactory; VarPtr fhat = varFactory.fluxVar("\\widehat{f}"); VarPtr u = varFactory.fieldVar("u"); VarPtr v = varFactory.testVar("v",HGRAD); BFPtr bf = Teuchos::rcp( new BF(varFactory) ); // initialize bilinear form //////////////////////////////////////////////////////////////////// // CREATE MESH //////////////////////////////////////////////////////////////////// // create a pointer to a new mesh: Teuchos::RCP<Mesh> mesh = Mesh::buildQuadMesh(meshPoints, horizontalCells, verticalCells, bf, H1Order, H1Order+pToAdd, useTriangles); mesh->setPartitionPolicy(Teuchos::rcp(new ZoltanMeshPartitionPolicy("HSFC"))); //////////////////////////////////////////////////////////////////// // INITIALIZE BACKGROUND FLOW FUNCTIONS //////////////////////////////////////////////////////////////////// BCPtr nullBC = Teuchos::rcp((BC*)NULL); RHSPtr nullRHS = Teuchos::rcp((RHS*)NULL); IPPtr nullIP = Teuchos::rcp((IP*)NULL); SolutionPtr backgroundFlow = Teuchos::rcp(new Solution(mesh, nullBC, nullRHS, nullIP) ); 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 ) ); //////////////////////////////////////////////////////////////////// // DEFINE BILINEAR FORM //////////////////////////////////////////////////////////////////// // v: // (sigma, grad v)_K - (sigma_hat_n, v)_dK - (u, beta dot grad v) + (u_hat * n dot beta, v)_dK bf->addTerm( -u, beta * v->grad()); bf->addTerm( fhat, v); // ==================== SET INITIAL GUESS ========================== mesh->registerSolution(backgroundFlow); FunctionPtr zero = Teuchos::rcp( new ConstantScalarFunction(0.0) ); FunctionPtr u0 = Teuchos::rcp( new U0 ); map<int, Teuchos::RCP<Function> > functionMap; functionMap[u->ID()] = u0; backgroundFlow->projectOntoMesh(functionMap); // ==================== END SET INITIAL GUESS ========================== //////////////////////////////////////////////////////////////////// // DEFINE INNER PRODUCT //////////////////////////////////////////////////////////////////// IPPtr ip = Teuchos::rcp( new IP ); ip->addTerm( v ); ip->addTerm( beta * v->grad() ); //////////////////////////////////////////////////////////////////// // DEFINE RHS //////////////////////////////////////////////////////////////////// Teuchos::RCP<RHSEasy> rhs = Teuchos::rcp( new RHSEasy ); FunctionPtr u_prev_squared_div2 = 0.5 * u_prev * u_prev; rhs->addTerm( (e1 * u_prev_squared_div2 + e2 * u_prev) * v->grad()); //////////////////////////////////////////////////////////////////// // DEFINE DIRICHLET BC //////////////////////////////////////////////////////////////////// Teuchos::RCP<BCEasy> inflowBC = Teuchos::rcp( new BCEasy ); // Create spatial filters SpatialFilterPtr bottomBoundary = Teuchos::rcp( new BottomBoundary ); SpatialFilterPtr leftBoundary = Teuchos::rcp( new LeftBoundary ); SpatialFilterPtr rightBoundary = Teuchos::rcp( new LeftBoundary ); // Create BCs FunctionPtr n = Teuchos::rcp( new UnitNormalFunction ); FunctionPtr u0_squared_div_2 = 0.5 * u0 * u0; SimpleFunction* u0Ptr = static_cast<SimpleFunction *>(u0.get()); double u0Left = u0Ptr->value(0,0); double u0Right = u0Ptr->value(1.0,0); FunctionPtr leftVal = Teuchos::rcp( new ConstantScalarFunction( -0.5*u0Left*u0Left ) ); FunctionPtr rightVal = Teuchos::rcp( new ConstantScalarFunction( 0.5*u0Right*u0Right ) ); inflowBC->addDirichlet(fhat, bottomBoundary, -u0 ); inflowBC->addDirichlet(fhat, leftBoundary, leftVal ); inflowBC->addDirichlet(fhat, rightBoundary, rightVal ); //////////////////////////////////////////////////////////////////// // CREATE SOLUTION OBJECT //////////////////////////////////////////////////////////////////// Teuchos::RCP<Solution> solution = Teuchos::rcp(new Solution(mesh, inflowBC, rhs, ip)); mesh->registerSolution(solution); if (enforceLocalConservation) { FunctionPtr zero = Teuchos::rcp( new ConstantScalarFunction(0.0) ); solution->lagrangeConstraints()->addConstraint(fhat == zero); } //////////////////////////////////////////////////////////////////// // DEFINE REFINEMENT STRATEGY //////////////////////////////////////////////////////////////////// Teuchos::RCP<RefinementStrategy> refinementStrategy; refinementStrategy = Teuchos::rcp(new RefinementStrategy(solution,energyThreshold)); //////////////////////////////////////////////////////////////////// // SOLVE //////////////////////////////////////////////////////////////////// for (int refIndex=0; refIndex<=numRefs; refIndex++) { double L2Update = 1e7; int iterCount = 0; while (L2Update > nonlinearRelativeEnergyTolerance && iterCount < maxNewtonIterations) { solution->solve(); L2Update = solution->L2NormOfSolutionGlobal(u->ID()); cout << "L2 Norm of Update = " << L2Update << endl; // backgroundFlow->clear(); backgroundFlow->addSolution(solution, newtonStepSize); iterCount++; } cout << endl; // check conservation VarPtr testOne = varFactory.testVar("1", CONSTANT_SCALAR); // Create a fake bilinear form for the testing BFPtr fakeBF = Teuchos::rcp( new BF(varFactory) ); // Define our mass flux FunctionPtr massFlux = Teuchos::rcp( new PreviousSolutionFunction(solution, fhat) ); LinearTermPtr massFluxTerm = massFlux * testOne; Teuchos::RCP<shards::CellTopology> quadTopoPtr = Teuchos::rcp(new shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() )); DofOrderingFactory dofOrderingFactory(fakeBF); int fakeTestOrder = H1Order; DofOrderingPtr testOrdering = dofOrderingFactory.testOrdering(fakeTestOrder, *quadTopoPtr); int testOneIndex = testOrdering->getDofIndex(testOne->ID(),0); vector< ElementTypePtr > elemTypes = mesh->elementTypes(); // global element types map<int, double> massFluxIntegral; // cellID -> integral double maxMassFluxIntegral = 0.0; double totalMassFlux = 0.0; double totalAbsMassFlux = 0.0; for (vector< ElementTypePtr >::iterator elemTypeIt = elemTypes.begin(); elemTypeIt != elemTypes.end(); elemTypeIt++) { ElementTypePtr elemType = *elemTypeIt; vector< ElementPtr > elems = mesh->elementsOfTypeGlobal(elemType); vector<int> cellIDs; for (int i=0; i<elems.size(); i++) { cellIDs.push_back(elems[i]->cellID()); } FieldContainer<double> physicalCellNodes = mesh->physicalCellNodesGlobal(elemType); BasisCachePtr basisCache = Teuchos::rcp( new BasisCache(elemType,mesh) ); basisCache->setPhysicalCellNodes(physicalCellNodes,cellIDs,true); // true: create side caches FieldContainer<double> cellMeasures = basisCache->getCellMeasures(); FieldContainer<double> fakeRHSIntegrals(elems.size(),testOrdering->totalDofs()); massFluxTerm->integrate(fakeRHSIntegrals,testOrdering,basisCache,true); // true: force side evaluation for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; // pick out the ones for testOne: massFluxIntegral[cellID] = fakeRHSIntegrals(i,testOneIndex); } // find the largest: for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; maxMassFluxIntegral = max(abs(massFluxIntegral[cellID]), maxMassFluxIntegral); } for (int i=0; i<elems.size(); i++) { int cellID = cellIDs[i]; maxMassFluxIntegral = max(abs(massFluxIntegral[cellID]), maxMassFluxIntegral); totalMassFlux += massFluxIntegral[cellID]; totalAbsMassFlux += abs( massFluxIntegral[cellID] ); } } if (rank==0) { cout << endl; cout << "largest mass flux: " << maxMassFluxIntegral << endl; cout << "total mass flux: " << totalMassFlux << endl; cout << "sum of mass flux absolute value: " << totalAbsMassFlux << endl; cout << endl; stringstream outfile; outfile << "burgers_" << refIndex; backgroundFlow->writeToVTK(outfile.str(), 5); } if (refIndex < numRefs) refinementStrategy->refine(rank==0); // print to console on rank 0 } return 0; }