libMesh
Functions
adjoints_ex6.C File Reference

Go to the source code of this file.

Functions

void write_output (EquationSystems &es, unsigned int a_step, std::string solution_type)
 
void set_system_parameters (PoissonSystem &system, FEMParameters &param)
 
std::unique_ptr< MeshRefinementbuild_mesh_refinement (MeshBase &mesh, FEMParameters &param)
 
std::unique_ptr< AdjointRefinementEstimatorbuild_adjoint_refinement_error_estimator (QoISet &qois)
 
int main (int argc, char **argv)
 

Function Documentation

◆ build_adjoint_refinement_error_estimator()

std::unique_ptr<AdjointRefinementEstimator> build_adjoint_refinement_error_estimator ( QoISet qois)

Definition at line 204 of file adjoints_ex6.C.

References libMesh::out.

Referenced by main().

205 {
206  libMesh::out << "Computing the error estimate using the Adjoint Refinement Error Estimator\n" << std::endl;
207 
208  auto adjoint_refinement_estimator = std::make_unique<AdjointRefinementEstimator>();
209 
210  adjoint_refinement_estimator->qoi_set() = qois;
211 
212  // We enrich the FE space for the dual problem by doing 2 uniform h refinements
213  adjoint_refinement_estimator->number_h_refinements = 2;
214 
215  return adjoint_refinement_estimator;
216 }
OStreamProxy out

◆ build_mesh_refinement()

std::unique_ptr<MeshRefinement> build_mesh_refinement ( MeshBase mesh,
FEMParameters param 
)

Definition at line 184 of file adjoints_ex6.C.

References FEMParameters::coarsen_fraction, FEMParameters::coarsen_threshold, FEMParameters::global_tolerance, mesh, FEMParameters::nelem_target, and FEMParameters::refine_fraction.

Referenced by main().

186 {
187  auto mesh_refinement = std::make_unique<MeshRefinement>(mesh);
188  mesh_refinement->coarsen_by_parents() = true;
189  mesh_refinement->absolute_global_tolerance() = param.global_tolerance;
190  mesh_refinement->nelem_target() = param.nelem_target;
191  mesh_refinement->refine_fraction() = param.refine_fraction;
192  mesh_refinement->coarsen_fraction() = param.coarsen_fraction;
193  mesh_refinement->coarsen_threshold() = param.coarsen_threshold;
194 
195  return mesh_refinement;
196 }
unsigned int nelem_target
Definition: femparameters.h:60
MeshBase & mesh
libMesh::Real refine_fraction
Definition: femparameters.h:62
libMesh::Real coarsen_threshold
Definition: femparameters.h:62
libMesh::Real global_tolerance
Definition: femparameters.h:61
libMesh::Real coarsen_fraction
Definition: femparameters.h:62

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 222 of file adjoints_ex6.C.

References libMesh::QoISet::add_indices(), libMesh::EquationSystems::add_system(), libMesh::DifferentiableSystem::adjoint_solve(), libMesh::MeshBase::all_second_order(), libMesh::DifferentiableQoI::assemble_qoi_sides, build_adjoint_refinement_error_estimator(), build_mesh_refinement(), libMesh::MeshTools::Generation::build_square(), libMesh::default_solver_package(), libMesh::err, libMesh::System::get_adjoint_solution(), libMesh::DifferentiableSystem::get_linear_solver(), PoissonSystem::get_QoI_value(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::EquationSystems::n_active_dofs(), libMesh::MeshBase::n_active_elem(), libMesh::out, PoissonSystem::postprocess(), libMesh::DifferentiableSystem::postprocess_sides, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::ParallelObject::processor_id(), libMesh::QUAD4, FEMParameters::read(), libMesh::EquationSystems::reinit(), libMesh::LinearSolver< T >::reuse_preconditioner(), libMesh::System::set_adjoint_already_solved(), set_system_parameters(), libMesh::QoISet::set_weight(), libMesh::System::solution, libMesh::FEMSystem::solve(), libMesh::NumericVector< T >::swap(), libMesh::MeshRefinement::uniformly_refine(), and write_output().

223 {
224  // Initialize libMesh.
225  LibMeshInit init (argc, argv);
226 
227  // This example requires a linear solver package.
228  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
229  "--enable-petsc, --enable-trilinos, or --enable-eigen");
230 
231  // Skip adaptive examples on a non-adaptive libMesh build
232 #ifndef LIBMESH_ENABLE_AMR
233  libmesh_example_requires(false, "--enable-amr");
234 #else
235 
236  // We use Dirichlet boundary conditions here
237 #ifndef LIBMESH_ENABLE_DIRICHLET
238  libmesh_example_requires(false, "--enable-dirichlet");
239 #endif
240 
241  libMesh::out << "Started " << argv[0] << std::endl;
242 
243  // Make sure the general input file exists, and parse it
244  {
245  std::ifstream i("general.in");
246  libmesh_error_msg_if(!i, '[' << init.comm().rank() << "] Can't find general.in; exiting early.");
247  }
248 
249  // Read in parameters from the input file
250  GetPot infile("general.in");
251 
252  // But allow the command line to override it.
253  infile.parse_command_line(argc, argv);
254 
255  FEMParameters param(init.comm());
256  param.read(infile);
257 
258  // Skip this default-2D example if libMesh was compiled as 1D-only.
259  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
260 
261  // Create a mesh, with dimension to be overridden later, distributed
262  // across the default MPI communicator.
263  Mesh mesh(init.comm());
264 
265  // And an object to refine it
266  std::unique_ptr<MeshRefinement> mesh_refinement =
267  build_mesh_refinement(mesh, param);
268 
269  // And an EquationSystems to run on it
270  EquationSystems equation_systems (mesh);
271 
272  libMesh::out << "Reading in and building the mesh" << std::endl;
273 
274  // Read in the mesh
276  (mesh, 2, 2,
277  0.0, 1.0,
278  0.0, 1.0,
279  QUAD4);
280 
281  // Make all the elements of the mesh second order so we can compute
282  // with a higher order basis
284 
285  // Create a mesh refinement object to do the initial uniform refinements
286  // on the coarse grid read in from lshaped.xda
287  MeshRefinement initial_uniform_refinements(mesh);
288  initial_uniform_refinements.uniformly_refine(param.coarserefinements);
289 
290  libMesh::out << "Building system" << std::endl;
291 
292  // Build the FEMSystem
293  PoissonSystem & system = equation_systems.add_system<PoissonSystem> ("PoissonSystem");
294 
295  // Set its parameters
296  set_system_parameters(system, param);
297 
298  libMesh::out << "Initializing systems" << std::endl;
299 
300  equation_systems.init ();
301 
302  // Print information about the mesh and system to the screen.
303  mesh.print_info();
304  equation_systems.print_info();
305 
306  // Get a pointer to the linear solver object to be able to reuse preconditioner
307  LinearSolver<Number> * linear_solver = system.get_linear_solver();
308 
309  {
310  // Adaptively solve the timestep
311  unsigned int a_step = 0;
312  for (; a_step != param.max_adaptivesteps; ++a_step)
313  {
314  // We can't adapt to both a tolerance and a
315  // target mesh size
316  if (param.global_tolerance != 0.)
317  libmesh_assert_equal_to (param.nelem_target, 0);
318  // If we aren't adapting to a tolerance we need a
319  // target mesh size
320  else
321  libmesh_assert_greater (param.nelem_target, 0);
322 
323  // Dont reuse preconditioners before the primal solve
324  linear_solver->reuse_preconditioner(false);
325 
326  // Solve the forward problem
327  system.solve();
328 
329  // Write out the computed primal solution
330  write_output(equation_systems, a_step, "primal");
331 
332  // Declare a QoISet object, we need this object to set weights for our QoI error contributions
333  QoISet qois;
334 
335  // Each index will correspond to a QoI
336  qois.add_indices({0});
337 
338  // Set weights for each index, these will weight the contribution of each QoI in the final error
339  // estimate to be used for flagging elements for refinement
340  qois.set_weight(0, 1.0);
341 
342  // Make sure we get the contributions to the adjoint RHS from the sides
343  system.assemble_qoi_sides = true;
344 
345  // We are about to solve the adjoint system, but before we do this we see the same preconditioner
346  // flag to reuse the preconditioner from the forward solver
347  linear_solver->reuse_preconditioner(param.reuse_preconditioner);
348 
349  // Solve the adjoint system. This takes the transpose of the stiffness matrix and then
350  // solves the resulting system
351  system.adjoint_solve();
352 
353  //Now that we have solved the adjoint, set the adjoint_already_solved boolean to true, so we dont solve unnecessarily in the error estimator
354  system.set_adjoint_already_solved(true);
355 
356  // Get a pointer to the primal solution vector
357  NumericVector<Number> & primal_solution = *system.solution;
358 
359  //Get a pointer to the solution vector of the adjoint problem for QoI 0
360  NumericVector<Number> & dual_solution_0 = system.get_adjoint_solution(0);
361 
362  //Swap the primal and dual solutions so we can write out the adjoint solution
363  primal_solution.swap(dual_solution_0);
364  write_output(equation_systems, a_step, "adjoint_0");
365 
366  //Swap back
367  primal_solution.swap(dual_solution_0);
368 
369  libMesh::out << "Adaptive step " << a_step << ", we have " << mesh.n_active_elem()
370  << " active elements and "
371  << equation_systems.n_active_dofs()
372  << " active dofs." << std::endl ;
373 
374  // Postprocess, compute the approximate QoIs and write them out to the console
375  libMesh::out << "Postprocessing: " << std::endl;
376  system.postprocess_sides = true;
377  system.postprocess();
378 
379  Number QoI_0_computed = system.get_QoI_value("computed", 0);
380  Number QoI_0_exact = system.get_QoI_value("exact", 0);
381 
382  libMesh::out << "The computed QoI 0 is " << std::setprecision(17)
383  << QoI_0_computed << std::endl;
384  libMesh::out << "The relative error in QoI 0 is " << std::setprecision(17)
385  << std::abs(QoI_0_computed - QoI_0_exact) << std::endl; // / std::abs(QoI_0_exact)
386 
387  // We will declare an error vector for passing to the adjoint refinement error estimator
388  ErrorVector QoI_elementwise_error;
389 
390  // Build an adjoint refinement error estimator object
391  std::unique_ptr<AdjointRefinementEstimator> adjoint_refinement_error_estimator =
393 
394  // Estimate the error in each element using the Adjoint Refinement estimator
395  adjoint_refinement_error_estimator->estimate_error(system, QoI_elementwise_error);
396 
397  // Print out the computed error estimate, note that we access the global error estimates
398  // using an accessor function, right now sum(QoI_elementwise_error) != global_QoI_error_estimate
399  libMesh::out << "The computed relative error in QoI 0 is " << std::setprecision(17)
400  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) << std::endl; // / std::abs(QoI_0_exact)
401 
402  // Also print out effectivity indices (estimated error/true error)
403  libMesh::out << "The effectivity index for the computed error in QoI 0 is " << std::setprecision(17)
404  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) /
405  std::abs(QoI_0_computed - QoI_0_exact) << std::endl;
406 
407  // For refinement purposes we need to sort by error
408  // *magnitudes*, but AdjointRefinement gives us signed errors.
409  if (!param.refine_uniformly)
410  for (std::size_t i=0; i<QoI_elementwise_error.size(); i++)
411  if (QoI_elementwise_error[i] != 0.)
412  QoI_elementwise_error[i] = std::abs(QoI_elementwise_error[i]);
413 
414  // We have to refine either based on reaching an error tolerance or
415  // a number of elements target, which should be verified above
416  // Otherwise we flag elements by error tolerance or nelem target
417 
418  // Uniform refinement
419  if (param.refine_uniformly)
420  {
421  mesh_refinement->uniformly_refine(1);
422  }
423  // Adaptively refine based on reaching an error tolerance
424  else if (param.global_tolerance >= 0. && param.nelem_target == 0.)
425  {
426  mesh_refinement->flag_elements_by_error_tolerance (QoI_elementwise_error);
427 
428  mesh_refinement->refine_and_coarsen_elements();
429  }
430  // Adaptively refine based on reaching a target number of elements
431  else
432  {
433  if (mesh.n_active_elem() >= param.nelem_target)
434  {
435  libMesh::out << "We reached the target number of elements.\n" << std::endl;
436  break;
437  }
438 
439  mesh_refinement->flag_elements_by_nelem_target (QoI_elementwise_error);
440 
441  mesh_refinement->refine_and_coarsen_elements();
442  }
443 
444  // Dont forget to reinit the system after each adaptive refinement!
445  equation_systems.reinit();
446 
447  libMesh::out << "Refined mesh to "
448  << mesh.n_active_elem()
449  << " active elements and "
450  << equation_systems.n_active_dofs()
451  << " active dofs." << std::endl;
452  }
453 
454  // On the last adaptive step, dont refine elements and check regressions via asserts
455  if (a_step == param.max_adaptivesteps)
456  {
457  linear_solver->reuse_preconditioner(false);
458  system.solve();
459 
460  write_output(equation_systems, a_step, "primal");
461 
462  NumericVector<Number> & primal_solution = *system.solution;
463 
464  QoISet qois;
465  std::vector<unsigned int> qoi_indices;
466 
467  qoi_indices.push_back(0);
468  qois.add_indices(qoi_indices);
469 
470  qois.set_weight(0, 1.0);
471 
472  system.assemble_qoi_sides = true;
473 
474  linear_solver->reuse_preconditioner(param.reuse_preconditioner);
475  system.adjoint_solve();
476  system.set_adjoint_already_solved(true);
477 
478  NumericVector<Number> & dual_solution_0 = system.get_adjoint_solution(0);
479 
480  primal_solution.swap(dual_solution_0);
481  write_output(equation_systems, a_step, "adjoint_0");
482 
483  primal_solution.swap(dual_solution_0);
484 
485  libMesh::out << "Adaptive step " << a_step << ", we have " << mesh.n_active_elem()
486  << " active elements and "
487  << equation_systems.n_active_dofs()
488  << " active dofs." << std::endl ;
489 
490  libMesh::out << "Postprocessing: " << std::endl;
491  system.postprocess_sides = true;
492  system.postprocess();
493 
494  Number QoI_0_computed = system.get_QoI_value("computed", 0);
495  Number QoI_0_exact = system.get_QoI_value("exact", 0);
496 
497  libMesh::out << "The computed QoI 0 is " << std::setprecision(17)
498  << QoI_0_computed << std::endl;
499  libMesh::out << "The relative error in QoI 0 is " << std::setprecision(17)
500  << std::abs(QoI_0_computed - QoI_0_exact) << std::endl; // / std::abs(QoI_0_exact)
501 
502 
503  // We will declare an error vector for passing to the adjoint
504  // refinement error estimator Right now, only the first entry
505  // of this vector will be filled (with the global QoI error
506  // estimate) Later, each entry of the vector will contain
507  // elementwise error that the user can sum to get the total
508  // error
509  ErrorVector QoI_elementwise_error;
510 
511  // Build an adjoint refinement error estimator object
512  std::unique_ptr<AdjointRefinementEstimator> adjoint_refinement_error_estimator =
514 
515  // Estimate the error in each element using the Adjoint Refinement estimator
516  adjoint_refinement_error_estimator->estimate_error(system, QoI_elementwise_error);
517 
518  // Print out the computed error estimate, note that we access
519  // the global error estimates using an accessor function,
520  // right now sum(QoI_elementwise_error) != global_QoI_error_estimate
521  libMesh::out << "The computed relative error in QoI 0 is " << std::setprecision(17)
522  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) << std::endl; // / std::abs(QoI_0_exact)
523 
524  // Also print out effectivity indices (estimated error/true error)
525  libMesh::out << "The effectivity index for the computed error in QoI 0 is " << std::setprecision(17)
526  << std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) /
527  std::abs(QoI_0_computed - QoI_0_exact) << std::endl;
528 
529  // Hard coded assert to ensure that the actual numbers we are getting are what they should be
530 
531  // The effectivity index isn't exactly reproducible at single precision
532  // libmesh_assert_less(std::abs(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact) - 0.84010976704434637), 1.e-5);
533  // libmesh_assert_less(std::abs(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact) - 0.48294428289950514), 1.e-5);
534 
535  // But the effectivity indices should always be sane
536  // libmesh_assert_less(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact), 2.5);
537  // libmesh_assert_greater(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(0)) / std::abs(QoI_0_computed - QoI_0_exact), .4);
538  // libmesh_assert_less(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact), 2.5);
539  // libmesh_assert_greater(std::abs(adjoint_refinement_error_estimator->get_global_QoI_error_estimate(1)) / std::abs(QoI_1_computed - QoI_1_exact), .4);
540 
541  // And the computed errors should still be low
542  // libmesh_assert_less(std::abs(QoI_0_computed - QoI_0_exact), 2e-4);
543  // libmesh_assert_less(std::abs(QoI_1_computed - QoI_1_exact), 2e-4);
544  }
545  }
546 
547  libMesh::err << '[' << mesh.processor_id()
548  << "] Completing output." << std::endl;
549 
550 #endif // #ifndef LIBMESH_ENABLE_AMR
551 
552  // All done.
553  return 0;
554 }
OStreamProxy err
This is the EquationSystems class.
void set_weight(std::size_t, Real)
Set the weight for this index.
Definition: qoi_set.h:232
virtual dof_id_type n_active_elem() const =0
void set_adjoint_already_solved(bool setting)
Setter for the adjoint_already_solved boolean.
Definition: system.h:415
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
virtual void postprocess(void)
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides...
Definition: poisson.C:182
MeshBase & mesh
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
bool postprocess_sides
If postprocess_sides is true (it is false by default), the postprocessing loop will loop over all sid...
Definition: diff_system.h:334
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
SolverPackage default_solver_package()
Definition: libmesh.C:1117
Implements (adaptive) mesh refinement algorithms for a MeshBase.
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1562
void write_output(EquationSystems &es, unsigned int a_step, std::string solution_type)
Definition: adjoints_ex6.C:98
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1593
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
void read(GetPot &input, const std::vector< std::string > *other_variable_names=nullptr)
std::unique_ptr< MeshRefinement > build_mesh_refinement(MeshBase &mesh, FEMParameters &param)
Definition: adjoints_ex6.C:184
virtual LinearSolver< Number > * get_linear_solver() const override
Definition: diff_system.C:163
virtual void reuse_preconditioner(bool)
Set the same_preconditioner flag, which indicates if we reuse the same preconditioner for subsequent ...
virtual void swap(NumericVector< T > &v)
Swaps the contents of this with v.
bool assemble_qoi_sides
If assemble_qoi_sides is true (it is false by default), the assembly loop for a quantity of interest ...
Definition: diff_qoi.h:99
OStreamProxy out
void add_indices(const std::vector< unsigned int > &indices)
Add this indices to the set to be calculated.
Definition: qoi_set.C:46
void all_second_order(const bool full_ordered=true)
Calls the range-based version of this function with a range consisting of all elements in the mesh...
Definition: mesh_base.C:1608
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
Definition: system.C:1245
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1076
void set_system_parameters(PoissonSystem &system, FEMParameters &param)
Definition: adjoints_ex6.C:122
std::unique_ptr< AdjointRefinementEstimator > build_adjoint_refinement_error_estimator(QoISet &qois)
Definition: adjoints_ex6.C:204
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
processor_id_type processor_id() const
Number & get_QoI_value(std::string type, unsigned int QoI_index)
Definition: poisson.h:29
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices=QoISet()) override
This function sets the _is_adjoint boolean member of TimeSolver to true and then calls the adjoint_so...
Definition: diff_system.C:150

◆ set_system_parameters()

void set_system_parameters ( PoissonSystem system,
FEMParameters param 
)

Definition at line 122 of file adjoints_ex6.C.

References PoissonSystem::analytic_jacobians(), FEMParameters::analytic_jacobians, FEMParameters::constrain_in_solver, PoissonSystem::fe_family(), FEMParameters::fe_family, PoissonSystem::fe_order(), FEMParameters::fe_order, FEMParameters::initial_linear_tolerance, FEMParameters::linear_tolerance_multiplier, FEMParameters::max_linear_iterations, FEMParameters::max_nonlinear_iterations, FEMParameters::min_step_length, FEMParameters::minimum_linear_tolerance, FEMParameters::print_jacobian_norms, libMesh::DifferentiableSystem::print_jacobian_norms, FEMParameters::print_jacobians, libMesh::DifferentiableSystem::print_jacobians, FEMParameters::print_residual_norms, libMesh::DifferentiableSystem::print_residual_norms, FEMParameters::print_residuals, libMesh::DifferentiableSystem::print_residuals, FEMParameters::print_solution_norms, libMesh::DifferentiableSystem::print_solution_norms, FEMParameters::print_solutions, libMesh::DifferentiableSystem::print_solutions, FEMParameters::relative_residual_tolerance, FEMParameters::relative_step_tolerance, FEMParameters::require_residual_reduction, libMesh::DifferentiableSystem::set_constrain_in_solver(), FEMParameters::solver_quiet, libMesh::DifferentiableSystem::time_solver, FEMParameters::use_petsc_snes, FEMParameters::verify_analytic_jacobians, and libMesh::FEMSystem::verify_analytic_jacobians.

Referenced by main().

123 {
124  // Use analytical jacobians?
125  system.analytic_jacobians() = param.analytic_jacobians;
126 
127  // Verify analytic jacobians against numerical ones?
129 
130  // Use the prescribed FE type
131  system.fe_family() = param.fe_family[0];
132  system.fe_order() = param.fe_order[0];
133 
134  // More desperate debugging options
136  system.print_solutions = param.print_solutions;
138  system.print_residuals = param.print_residuals;
140  system.print_jacobians = param.print_jacobians;
141 
142  // No transient time solver
143  system.time_solver = std::make_unique<SteadySolver>(system);
144 
145  // Nonlinear solver options
146  if (param.use_petsc_snes)
147  {
148 #ifdef LIBMESH_HAVE_PETSC
149  system.time_solver->diff_solver() = std::make_unique<PetscDiffSolver>(system);
150 #else
151  libmesh_error_msg("This example requires libMesh to be compiled with PETSc support.");
152 #endif
153  }
154  else
155  {
156  system.time_solver->diff_solver() = std::make_unique<NewtonSolver>(system);
157  auto solver = cast_ptr<NewtonSolver*>(system.time_solver->diff_solver().get());
158 
159  solver->quiet = param.solver_quiet;
160  solver->max_nonlinear_iterations = param.max_nonlinear_iterations;
161  solver->minsteplength = param.min_step_length;
162  solver->relative_step_tolerance = param.relative_step_tolerance;
163  solver->relative_residual_tolerance = param.relative_residual_tolerance;
164  solver->require_residual_reduction = param.require_residual_reduction;
165  solver->linear_tolerance_multiplier = param.linear_tolerance_multiplier;
166  if (system.time_solver->reduce_deltat_on_diffsolver_failure)
167  {
168  solver->continue_after_max_iterations = true;
169  solver->continue_after_backtrack_failure = true;
170  }
172 
173  // And the linear solver options
174  solver->max_linear_iterations = param.max_linear_iterations;
175  solver->initial_linear_tolerance = param.initial_linear_tolerance;
176  solver->minimum_linear_tolerance = param.minimum_linear_tolerance;
177  }
178 }
bool print_solution_norms
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
Definition: fem_system.h:215
bool analytic_jacobians
virtual void set_constrain_in_solver(bool enable)
set_constrain_in_solver to false to apply constraints only via residual terms in the systems to be so...
Definition: diff_system.C:424
bool print_jacobian_norms
Set print_jacobian_norms to true to print |J| whenever it is assembled.
Definition: diff_system.h:361
bool require_residual_reduction
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we&#39;re going to use.
Definition: diff_system.h:253
bool print_residual_norms
libMesh::Real relative_step_tolerance
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
Definition: diff_system.h:366
std::vector< unsigned int > fe_order
bool print_solution_norms
Set print_residual_norms to true to print |U| whenever it is used in an assembly() call...
Definition: diff_system.h:340
unsigned int max_linear_iterations
bool constrain_in_solver
bool print_residual_norms
Set print_residual_norms to true to print |F| whenever it is assembled.
Definition: diff_system.h:351
double minimum_linear_tolerance
std::vector< std::string > fe_family
libMesh::Real verify_analytic_jacobians
bool print_residuals
Set print_residuals to true to print F whenever it is assembled.
Definition: diff_system.h:356
double initial_linear_tolerance
std::string & fe_family()
Definition: poisson.h:21
bool print_solutions
Set print_solutions to true to print U whenever it is used in an assembly() call. ...
Definition: diff_system.h:346
libMesh::Real relative_residual_tolerance
unsigned int & fe_order()
Definition: poisson.h:22
bool print_jacobian_norms
libMesh::Real min_step_length
unsigned int max_nonlinear_iterations
double linear_tolerance_multiplier
bool & analytic_jacobians()
Definition: poisson.h:23

◆ write_output()

void write_output ( EquationSystems es,
unsigned int  a_step,
std::string  solution_type 
)

Definition at line 98 of file adjoints_ex6.C.

References libMesh::EquationSystems::get_mesh(), libMesh::libmesh_ignore(), mesh, and libMesh::MeshOutput< MT >::write_equation_systems().

Referenced by main().

101 {
102  // Avoid unused variable warnings in the --disable-gmv configuration
103  libmesh_ignore(es, a_step, solution_type);
104 
105 #ifdef LIBMESH_HAVE_GMV
106  MeshBase & mesh = es.get_mesh();
107 
108  std::ostringstream file_name_gmv;
109  file_name_gmv << solution_type
110  << ".out.gmv."
111  << std::setw(2)
112  << std::setfill('0')
113  << std::right
114  << a_step;
115 
117  (file_name_gmv.str(), es);
118 #endif
119 }
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
MeshBase & mesh
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:46
This is the MeshBase class.
Definition: mesh_base.h:75
void libmesh_ignore(const Args &...)
const MeshBase & get_mesh() const