void LinearPDEConstrainedObj::solveState(const Vector<double>& x) const { Tabs tab(0); PLAYA_MSG2(verb(), tab << "solving state"); PLAYA_MSG3(verb(), tab << "|x|=" << x.norm2()); PLAYA_MSG5(verb(), tab << "x=" << endl << tab << x.norm2()); setDiscreteFunctionVector(designVarVal(), x); /* solve the state equations in order */ for (int i=0; i<stateProbs_.size(); i++) { SolverState<double> status = stateProbs_[i].solve(solvers_[i], stateVarVals(i)); TEUCHOS_TEST_FOR_EXCEPTION(status.finalState() != SolveConverged, std::runtime_error, "state equation could not be solved: status=" << status.stateDescription()); } PLAYA_MSG2(verb(), tab << "done state solve"); /* do postprocessing */ statePostprocCallback(); }
void LinearPDEConstrainedObj ::solveStateAndAdjoint(const Vector<double>& x) const { Tabs tab(0); PLAYA_MSG2(verb(), tab << "solving state and adjoint"); PLAYA_MSG3(verb(), tab << "|x|=" << x.norm2()); PLAYA_MSG5(verb(), tab << "x=" << endl << tab << x.norm2()); Tabs tab1; setDiscreteFunctionVector(designVarVal(), x); PLAYA_MSG3(verb(), tab1 << "solving state eqns"); /* solve the state equations in order */ for (int i=0; i<stateProbs_.size(); i++) { SolverState<double> status = stateProbs_[i].solve(solvers_[i], stateVarVals(i)); /* if the solve failed, write out the design var and known state * variables */ if (status.finalState() != SolveConverged) { FieldWriter w = new VTKWriter("badSolve"); w.addMesh(Lagrangian().mesh()); w.addField("designVar", new ExprFieldWrapper(designVarVal())); for (int j=0; j<i; j++) { Expr tmp = stateVarVals(j).flatten(); for (int k=0; k<tmp.size(); k++) { w.addField("stateVar-"+Teuchos::toString(j)+"-"+Teuchos::toString(k), new ExprFieldWrapper(tmp[k])); } } w.write(); } TEUCHOS_TEST_FOR_EXCEPTION(status.finalState() != SolveConverged, std::runtime_error, "state equation " << i << " could not be solved: status=" << status.stateDescription()); } PLAYA_MSG3(verb(), tab1 << "done solving state eqns"); /* do postprocessing */ statePostprocCallback(); PLAYA_MSG3(verb(), tab1 << "solving adjoint eqns"); /* solve the adjoint equations in reverse order */ for (int i=adjointProbs_.size()-1; i>=0; i--) { SolverState<double> status = adjointProbs_[i].solve(solvers_[i], adjointVarVals(i)); /* if the solve failed, write out the design var and known state * and adjoint variables */ if (status.finalState() != SolveConverged) { FieldWriter w = new VTKWriter("badSolve"); w.addMesh(Lagrangian().mesh()); w.addField("designVar", new ExprFieldWrapper(designVarVal())); for (int j=0; j<stateProbs_.size(); j++) { Expr tmp = stateVarVals(j).flatten(); for (int k=0; k<tmp.size(); k++) { w.addField("stateVar-"+Teuchos::toString(j)+"-"+Teuchos::toString(k), new ExprFieldWrapper(tmp[k])); } } for (int j=adjointProbs_.size()-1; j>i; j--) { Expr tmp = adjointVarVals(j).flatten(); for (int k=0; k<tmp.size(); k++) { w.addField("adjointVar-"+Teuchos::toString(j)+"-"+Teuchos::toString(k), new ExprFieldWrapper(tmp[k])); } } w.write(); } TEUCHOS_TEST_FOR_EXCEPTION(status.finalState() != SolveConverged, std::runtime_error, "adjoint equation " << i << " could not be solved: status=" << status.stateDescription()); } PLAYA_MSG3(verb(), tab1 << "done solving adjoint eqns"); PLAYA_MSG2(verb(), tab1 << "done solving state and adjoint eqns"); }
int main(int argc, char** argv) { try { int nx = 32; double convTol = 1.0e-8; double lambda = 0.5; Sundance::setOption("nx", nx, "Number of elements"); Sundance::setOption("tol", convTol, "Convergence tolerance"); Sundance::setOption("lambda", lambda, "Lambda (parameter in Bratu's equation)"); Sundance::init(&argc, &argv); Out::root() << "Bratu problem (lambda=" << lambda << ")" << endl; Out::root() << "Fixed-point iteration" << endl << endl; VectorType<double> vecType = new EpetraVectorType(); MeshType meshType = new BasicSimplicialMeshType(); MeshSource mesher = new PartitionedLineMesher(0.0, 1.0, nx, meshType); Mesh mesh = mesher.getMesh(); CellFilter interior = new MaximalCellFilter(); CellFilter sides = new DimensionalCellFilter(mesh.spatialDim()-1); CellFilter left = sides.subset(new CoordinateValueCellPredicate(0, 0.0)); CellFilter right = sides.subset(new CoordinateValueCellPredicate(0, 1.0)); BasisFamily basis = new Lagrange(1); Expr u = new UnknownFunction(basis, "u"); Expr v = new TestFunction(basis, "v"); Expr grad = gradient(1); Expr x = new CoordExpr(0); const double pi = 4.0*atan(1.0); Expr uExact = sin(pi*x); Expr R = pi*pi*uExact - lambda*exp(uExact); QuadratureFamily quad4 = new GaussianQuadrature(4); QuadratureFamily quad2 = new GaussianQuadrature(2); DiscreteSpace discSpace(mesh, basis, vecType); Expr uPrev = new DiscreteFunction(discSpace, 0.5); Expr uCur = copyDiscreteFunction(uPrev); Expr eqn = Integral(interior, (grad*u)*(grad*v) - v*lambda*exp(uPrev) - v*R, quad4); Expr h = new CellDiameterExpr(); Expr bc = EssentialBC(left+right, v*u/h, quad4); LinearProblem prob(mesh, eqn, bc, v, u, vecType); Expr normSqExpr = Integral(interior, pow(u-uPrev, 2.0), quad2); Functional normSqFunc(mesh, normSqExpr, vecType); FunctionalEvaluator normSqEval = normSqFunc.evaluator(u, uCur); LinearSolver<double> linSolver = LinearSolverBuilder::createSolver("amesos.xml"); Out::root() << "Fixed-point iteration" << endl; int maxIters = 20; Expr soln ; bool converged = false; for (int i=0; i<maxIters; i++) { /* solve for the next u */ prob.solve(linSolver, uCur); /* evaluate the norm of (uCur-uPrev) using * the FunctionalEvaluator defined above */ double deltaU = sqrt(normSqEval.evaluate()); Out::root() << "Iter=" << setw(3) << i << " ||Delta u||=" << setw(20) << deltaU << endl; /* check for convergence */ if (deltaU < convTol) { soln = uCur; converged = true; break; } /* get the vector from the current discrete function */ Vector<double> uVec = getDiscreteFunctionVector(uCur); /* copy the vector into the previous discrete function */ setDiscreteFunctionVector(uPrev, uVec); } TEUCHOS_TEST_FOR_EXCEPTION(!converged, std::runtime_error, "Fixed point iteration did not converge after " << maxIters << " iterations"); FieldWriter writer = new DSVWriter("FixedPointBratu.dat"); writer.addMesh(mesh); writer.addField("soln", new ExprFieldWrapper(soln[0])); writer.write(); Out::root() << "Converged!" << endl << endl; double L2Err = L2Norm(mesh, interior, soln-uExact, quad4); Out::root() << "L2 Norm of error: " << L2Err << endl; Sundance::passFailTest(L2Err, 1.5/((double) nx*nx)); } catch(exception& e) { Sundance::handleException(e); } Sundance::finalize(); }