PETSc version 3.17.4
DMAdaptorAdapt
Creates a new DM that is adapted to the problem
Synopsis
#include "petscdmadaptor.h"
PetscErrorCode DMAdaptorAdapt(DMAdaptor adaptor, Vec x, DMAdaptationStrategy strategy, DM *adm, Vec *ax)
Not collective
Input Parameters
| adaptor | - The DMAdaptor object
|
| x | - The global approximate solution
|
| strategy | - The adaptation strategy
|
Output Parameters
| adm | - The adapted DM
|
| ax | - The adapted solution
|
Options database keys
| -snes_adapt <strategy> | - initial, sequential, multigrid
|
| -adapt_gradient_view | - View the Clement interpolant of the solution gradient
|
| -adapt_hessian_view | - View the Clement interpolant of the solution Hessian
|
| -adapt_metric_view | - View the metric tensor for adaptive mesh refinement
|
Note: The available adaptation strategies are
1) Adapt the initial mesh until a quality metric, e.g., a priori error bound, is satisfied
2) Solve the problem on a series of adapted meshes until a quality metric, e.g. a posteriori error bound, is satisfied
3) Solve the problem on a hierarchy of adapted meshes generated to satisfy a quality metric using multigrid
See Also
DMAdaptorSetSolver(), DMAdaptorCreate(), DMAdaptorAdapt()
Level
intermediate
Location
src/snes/utils/dmadapt.c
Index of all SNES routines
Table of Contents for all manual pages
Index of all manual pages