Commit 8a9f959f authored by Müller, Felix's avatar Müller, Felix

Add documentation for initfile parameters

parent ebf66820
Pipeline #4923 passed with stage
in 56 minutes and 31 seconds
# Parameter file {: #group-initfile }
## Summary
Many functionalities of an `AMDiS` program can be controlled using a parameter file. This file is passed as first argument to the `AMDiS` executable. If we want to run the `ellipt.2d` example using the parameter file `init/ellipt.dat.2d` we can do so using the following console command in the examples directory:
```
./ellipt.2d init/ellipt.dat.2d [other arguments]
```
We will give a list of possible parameters at the [end of this document](#list-of-parameters-and-values).
## Structure of a parameter file
Parameters in the parameter file have keys in the form of strings separated by `->`. This allows grouping of parameters that are related to the same things. High-level classes like `ProblemStat` have a string constructor argument where a name can be provided that relates to an initfile key.
```c++
examples/ellipt.cc:
using Param = LagrangeBasis<Grid, 2>;
ProblemStat<Param> prob("ellipt"); // make a ProblemStat with the name 'ellipt'
```
```
examples/init/ellipt.dat.2d:
ellipt->mesh: elliptMesh % name of the mesh created by 'ellipt'
elliptMesh->num cells: 8 8 % size of the mesh 'elliptMesh'
ellipt->solver: pcg % solver type used by 'ellipt'
ellipt->solver->max iteration: 10000 % solver parameter related to 'ellipt->solver'
```
As we saw above, parameters are separated from the keys with the `:` character and can be strings and numbers. Anything following a `%` character is considered a comment.
## Reading from a parameter file
In the section above we learned that `AMDiS` classes can be given names that can then be used to specify parameters in the parameter file. In addition to the options `AMDiS` classes already provide we can add user-defined parameters. For this we simply have to add a `key: value` pair to the parameter file as we would do with the built-in options and read the value somewhere in our program using the function `Parameters::get<T>`:
```c++
// (1) Overwrites 'value' with the value related to 'key' if found
template <class T>
static void get(std::string const& key, T& value);
// (2) Returns a std::optional containing the value related to 'key' if found
template <class T>
static std::optional<T> get(std::string const& key);
```
We will now give an example of reading an `int` parameter called `myParameter` from an initfile that defaults to `3` if no parameter file containing the key is provided.
```
// example parameter file entry
myParameter: 42
```
```c++
// using (1)
int i = 3;
Parameters::get("myParameter", i);
// using (2)
int i = Parameters::get<int>("myParameter").value_or(3);
```
## Writing to a parameter file
It is possible to change to value of a parameter given in a parameter file using the function
```c++
template <class T>
static void set(std::string const& key, T const& value);
```
This can be useful if we want to override certain values within our program.
## List of parameters and values
In the following section `<prefix>` denotes a placeholder for the name given to the class instance. `<treepath>` is a (possibly empty) comma-separated list of integers specifying a [`TreePath`](../GlobalBasis#globalbasis-using-treepath).
### ProblemStat
key | type | values | default value
------------------------------|--------|--------------------------------|-------------
`<prefix>->mesh` | string | Prefix used by the mesh | `mesh`
`<prefix>->restore->grid` | string | Filename of the grid file | Not optional<sup>1</sup>
`<prefix>->restore->solution` | string | Filename of the solution file | Not optional<sup>1</sup>
`<prefix>->symmetry` | string | Symmetry of the stiffness matrix, one of:<br> `spd`,<br> `symmetric`,<br>`hermitian`,<br>`structurally_symmetric`,<br>`unknown` | `unknown`
`<prefix>->solver` | string | Solver type, see [Solvers and Preconditioners](#solvers-and-preconditioners) | `default`
`<prefix>->marker[<treepath>]`| string | Prefix for marker parameters, see [Markers](#markers) | None<sup>2</sup>
`<prefix>->output[<treepath>]`| string | Prefix for filewriter parameters, see [Filewriters](#filewriters)| None<sup>3</sup>
`<prefix>->output` | string | as `<prefix>->output[]` | None<sup>3</sup>
<sup>1</sup>Only required when `ProblemStat::restore` is called<br>
<sup>2</sup>No marker will be set up if the `->strategy` subkey is not specified<br>
<sup>3</sup>No filewriter will be set up if the `->format` subkey is not specified<br>
### Markers
key | type | values | default value
-------------------------------------|--------|----------------------------------|-------------
`<prefix>->strategy` | string | Marking strategy, one of:<br>`0` (no marking),<br>`1` (global refinement),<br>`2` (maximum strategy),<br>`3` (equidistribution strategy),<br>`4` (guaranteed error reduction strategy) | Not optional
`<prefix>->info` | int | Verbosity | `0`
`<prefix>->max refinement level` | int | Maximum refinement level | No maximum
`<prefix>->min refinement level` | int | Minimum refinement level | `0`
`<prefix>->p`<sup>123</sup> | double | Power in estimator norm | `2.0`
`<prefix>->MSGamma`<sup>1</sup> | double | Marking parameter | `0.5`
`<prefix>->MSGammaC`<sup>1</sup> | double | Marking parameter | `0.1`
`<prefix>->ESTheta`<sup>2</sup> | double | Marking parameter | `0.9`
`<prefix>->EsThetaC`<sup>2</sup> | double | Marking parameter | `0.2`
`<prefix>->GERSThetaStar`<sup>3</sup>| double | Marking parameter | `0.6`
`<prefix>->GERSNu`<sup>3</sup> | double | Marking parameter | `0.1`
`<prefix>->GERSThetaC`<sup>3</sup> | double | Marking parameter | `0.1`
<sup>1</sup>Only for maximum strategy<br>
<sup>2</sup>Only for equidistribution strategy<br>
<sup>3</sup>Only for guaranteed error reduction strategy<br>
### Solvers and Preconditioners
#### Parameters used by the SolverInfo class
key | type | values | default value
-------------------------------------|--------|----------------------------------|-------------
`<prefix>->info` | int | Solver verbosity level | `0`
`<prefix>->break if tolerance not reached`|bool| If `true` throw an error if tolerance could not be reached| `false`
#### PETSc
key | type | values | default value
--------------------------------|---------|----------------------------------|-------------
`<prefix>->info` | int | Verbosity, one of:<br>`0` (no convergence information),<br>`1` (minimal output, print residual every 10th iteration),<br>`2` (full convergence output) | `0`
`<prefix>->ksp` | string | Name of a PETSc Krylov method<sup>a</sup>| `default`
`<prefix>` | string | alternative to `<prefix>->ksp` when not set | as above
`<prefix>->max it` |PetscInt | Verbosity | `PETSC_DEFAULT`
`<prefix>->rtol` |PetscReal| Relative convergence tolerance<sup>b</sup>| `PETSC_DEFAULT`
`<prefix>->atol` |PetscReal| Absolute convergence tolerance<sup>b</sup>| `PETSC_DEFAULT`
`<prefix>->divtol` |PetscReal| Divergence tolerance<sup>b</sup> | `PETSC_DEFAULT`
`<prefix>->prefix` | string | The options prefix for the solver<sup>c</sup> | No prefix set
`<prefix>->mat solver` | string | Sets the software that is used to perform the factorization <sup>d</sup>| `umfpack` (sequential matrix) or `mumps` (otherwise)
`<prefix>->pc` | string | The preconditioner type<sup>e</sup>| `default`
`<prefix>->pc->mat solver`<sup>1</sup>|string| Solver used by the LU preconditioner<sup>d</sup>| None
`<prefix>->pc->prefix` | string | The options prefix for the preconditioner<sup>f</sup>|No prefix set
`<prefix>->scale`<sup>2</sup> |PetscReal| Damping factor used by richardson solver | `1.0`
`<prefix>->restart`<sup>3</sup> |PetscInt | Number of iterations until restart | `30`
`<prefix>->gramschmidt`<sup>3</sup>|string| Orthogonalization routine used by gmres solvers<sup>g</sup>: `classical` or `modified` | Use PETSc default
`<prefix>->restart`<sup>3</sup> | string | Type of iterative refinement to use in the classical Gram Schmidt orthogonalization<sup>h</sup>, one of:<br>`never` (never refine),<br>`ifneeded` (refine if needed),<br>`always` (always refine)| Use PETSc default
`<prefix>->ell`<sup>4</sup> |PetscReal| Number of search directions in BiCGStab_ell| Use PETSc default
`<prefix>->restart`<sup>5</sup> |PetscInt | Number of iterations until restart | Use PETSc default
<sup>1</sup>Only for `lu` preconditioner<br>
<sup>2</sup>Only for richardson solver<br>
<sup>3</sup>Only for gmres solver types<br>
<sup>4</sup>Only for BiCGStab_ell solver<br>
<sup>5</sup>Only for GCR solver<br>
<sup>a</sup>See [PETSc documentation for Krylov methods](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPType.html)<br>
<sup>b</sup>See [PETSc documentation for Krylov method tolerances](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetTolerances.html)<br>
<sup>c</sup>See [PETSc documentation for Krylov method options](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetOptionsPrefix.html)<br>
<sup>d</sup>See [PETSc documentation for preconditioners](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCFactorSetMatSolverType.html)<br>
<sup>e</sup>See [PETSc documentation for preconditioners](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCType.html)<br>
<sup>f</sup>See [PETSc documentation for preconditioners](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCSetOptionsPrefix.html)<br>
<sup>g</sup>See [PETSc documentation for GMRes solver](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGMRESSetOrthogonalization.html)<br>
<sup>h</sup>See [PETSc documentation for GMRes solver](https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPGMRESSetCGSRefinementType.html)<br>
#### EIGEN
key | type | values | default value
-------------------------------------|----------|-------------------------------------|----------
`<prefix>` | string | Solver type, one of:<br>`cg` (conjugate gradient solver),<br>`bicgstab` (bi-conjugate gradient stabilized solver),<br>`bcgs` (same as `bicgstab`),<br>`minres` (minimal residual algorithm for symmetric matrices),<br>`gmres` (generalized minimal residual algorithm),<br>`dgmres` (restarted gmres with deflation),<br>`default` (same as `gmres`),<br>`umfpack`<sup>1</sup> (sparse LU factorization based on UmfPack),<br>`superlu`<sup>2</sup> (sparse LU factorization based on SuperLU library),<br>`direct` (same as `umfpack` or `superlu`) | `default`
`<prefix>->reuse pattern` | bool | Reuse matrix pattern if `true` | `false`
`<prefix>->relative tolerance` |RealScalar| Tolerance threshold used by the stopping criteria| floating-point epsilon
`<prefix>->max iteration` | int | Maximum number of iterations | `500`
`<prefix>->restart`<sup>34</sup> | int | Number of iterations until restart | `30`
`<prefix>->num eigenvalues`<sup>4</sup>| int | Number of eigenvalues to deflate at each restart| `0`
`<prefix>->max num eigenvalues`<sup>4</sup>|int | Maximum size of the deflation subspace| `5`
`<prefix>->precon->name`<sup>5</sup> | string | Preconditioner, one of:<br>`no` (no preconditioner),<br>`diag` (diagonal preconditioner),<br>`jacobi` (diagonal preconditioner),<br>`ilu` (incomplete LU factorization),<br>`ic` (incomplete Choleski factorization) | `no`
`<prefix>->precon->ordering`<sup>6</sup>|string | Ordering used in Choleski factorization: `amd` or `natural` | `amd`
`<prefix>->precon->drop tolerance`<sup>7</sup>|double| Drop any element whose magnitude is less than this tolerance| `1e-15l`
`<prefix>->precon->fill factor`<sup>7</sup>|int | This is used to compute the number of largest elements to keep on each row| `10`
<sup>1</sup>Requires UmfPack<br>
<sup>2</sup>Requires SuperLU<br>
<sup>3</sup>Only for `gmres` solver<br>
<sup>4</sup>Only for `dgmres` solver<br>
<sup>5</sup>Does not apply to `umfpack` or `superlu` solvers<br>
<sup>6</sup>Only for incomplete Choleski factorization<br>
<sup>7</sup>Only for `ilu` solver<br>
#### MTL
key | type | values | default value
-------------------------------------|--------|----------------------------------|-------------
`<prefix>` | string | Solver type, one of:<br>`cg` (conjugate gradient method),<br>`cgs` (stabilized conjugate gradient method),<br>`bicg` (BiCG method),<br>`bcgs` (stabilized BiCG(1) method),<br>`bicgstab` (as `bcgs`),<br>`bcgs2` (stabilized BiCG(2) method),<br>`bicgstab2` (as `bcgs2`),<br>`qmr` (quasi-minimal residual method),<br>`tfqmr` (transposed-free quasi-minimal residual method),<br>`bcgsl` (stabilized BiCG(ell) method),<br>`bicgstab_ell` (as `bcgsl`),<br>`gmres` (generalized minimal residual method),<br>`fgmres` (flexible gmres),<br>`minres` (minimal residual method),<br>`idrs` (induced dimension reduction method),<br>`gcr` (generalized conjugate residual method),<br>`preonly` (just apply a preconditioner),<br>`default` (as `gmres`),<br>`umfpack`<sup>1</sup> (external UMFPACK solver),<br>`direct` (as `umfpack`) | `default`
`<prefix>->absolute tolerance`<sup>2</sup>| type of matrix values | Absolute solver tolerance: ∥b-Ax∥ < `absolute tolerance` | `0`
`<prefix>->relative tolerance`<sup>2</sup>| type of matrix values | Relative solver tolerance: ∥b-Ax∥ / ∥b-Ax<sub>0</sub>∥ < `relative tolerance` | `1.e-6`
`<prefix>->max iteration`<sup>2</sup>| int | Maximal number of iterations | `1000`
`<prefix>->print cycle`<sup>2</sup> | int | Print solver residuals every `print cycle` iterations | `100`
`<prefix>->store symbolic`<sup>3</sup>| bool | | `false`
`<prefix>->symmetric strategy`<sup>3</sup>|int| | `0`
`<prefix>->alloc init`<sup>3</sup> | double | | `0.7`
`<prefix>->ell`<sup>4</sup> | int | Parameter of stabilized BiCG(`ell`)| `3`
`<prefix>->restart`<sup>5</sup> | int | Maximal number of orthogonalized vectors| `30`
`<prefix>->orthogonalization`<sup>6</sup>|int | `1` (gram schmidt) or `2` (householder)| `1`
`<prefix>->s`<sup>7</sup> | int | Parameter of IDR(`s`) | `30`
`<prefix>->precon`<sup>2</sup> | string | Left preconditioner type, one of:<br>`diag` (diagonal preconditioner),<br>`jacobi` (as `diag`),<br>`masslumping` (mass lumping preconditioner),<br>`ilu` (incomplete LU preconditioner),<br>`ilu0` (as `ilu`),<br>`no` (no preconditioner),<br>`identity` (as `no`),<br>`solver` (use a solver as preconditioner),<br>`default` (as `no`),<br>`hypre`<sup>8</sup> (Hypre preconditioner) | `default`
`<prefix>->precon->solver`<sup>9</sup>| int | Solver type used for the preconditioner, same options as for `<prefix>`| `default`
<sup>1</sup>Requires UmfPack<br>
<sup>2</sup>Does not apply to `umfpack` solver<br>
<sup>3</sup>Only `umfpack` solver<br>
<sup>4</sup>Only `bcgsl` solver<br>
<sup>5</sup>Only `gmres`, `fgmres` and `gcr` solvers<br>
<sup>6</sup>Only `gmres` and `fgmres` solvers<br>
<sup>7</sup>Only `idrs` solver<br>
<sup>8</sup>Requires Hypre, see `amdis/linearalgebra/mtl/HyprePrecon.hpp` for subkeys and possible values<br>
<sup>9</sup>Only `solver` preconditioner, subkeys use the same values as above e.g. `<prefix>->precon->solver->max iteration: 500`<br>
#### ISTL
key | type | values | default value
-------------------------------------|--------|----------------------------------|-------------
`<prefix>` | string | Solver type, one of:<br>`cg` (conjugate gradient method),<br>`pcg` (generalized preconditioned conjugate gradient solver),<br>`fcg`<sup>1</sup> (accelerated flexible conjugate gradient method),<br>`cfcg`<sup>1</sup> (complete flexible conjugate gradient method),<br>`bcgs` (stabilized bi-conjugate gradient method),<br>`bicgstab` (as `bcgs`),<br>`default` (as `bcgs`),<br>`minres` (minimal residul method),<br>`gmres` (generalized minimal residual method),<br>`fgmres`<sup>1</sup> (flexible generalized minimal residual (FGMRes) method),<br>`umfpack`<sup>2</sup> (external UMFPACK solver),<br>`ldl`<sup>2</sup> (external LDL solver),<br>`spqr`<sup>2</sup> (external SQPR solver),<br>`cholmod`<sup>2</sup> (external Cholmod solver),<br>`superlu`<sup>2</sup> (external SuperLU solver),<br>`direct` (as `umfpack` or `superlu`) | `default`
`<prefix>->info` | int | Verbosity | `0`
`<prefix>->max iteration`<sup>3</sup>| int | Maximal number of solver iterations| `500`
`<prefix>->relative tolerance`<sup>3</sup>|type of matrix values| Relative break tolerance| `1.e-6`
`<prefix>->restart`<sup>4</sup> | int | Restart parameter | `30`
`<prefix>->reuse vector`<sup>5</sup> | bool | Reuse vectors in subsequent calls to apply if `true`| `true`
`<prefix>->category` | string | `sequential` (sequential solver),<br>`overlapping` (overlapping parallel solver),<br>`nonoverlapping` (nonoverlapping parallel solver),<br>`default` (chooses depending on grid and process count)| None<sup>9</sup>
`<prefix>->precon`<sup>3</sup> | string | Preconditioner type, one of:<br>`diag` (diagonal preconditioner),<br>`jacobi` (as `diag`),<br>`gs` (Gauss-Seidel preconditioner),<br> `gauss_seidel` (as `gs`),<br>`sor` (Successive Overrelaxation methods),<br>`ssor` (Symmetric Successive Overrelaxation methods),<br>`pssor` (A parallel SSOR preconditioner (requires overlap)),<br>`richardson` (Richardson methods),<br>`default` (as `richardson`),<br>`solver` (Turns a solver into a preconditioner),<br>`bjacobi` (Block-Jacobi methods),<br>`ilu` (Incomplete LU factorization),<br>`ilu0` (as `ilu`),<br>`ildl` (Incomplete LDL factorization),<br>`amg`<sup>6</sup> (Algebraic multigrid method),<br>`fastamg`<sup>6</sup> (Algebraic multigrid method),<br>`kamg`<sup>6</sup> (Algebraic multigrid method) | `default`
`<prefix>->precon->relaxation`<sup>3</sup>|double| Dumping/relaxation factor | `1.0`
`<prefix>->precon->iterations`<sup>3</sup>|int| Number of iterations | `1`
`<prefix>->precon->solver`<sup>7</sup>|string | Solver type used for the preconditioner, same options as for `<prefix>`| `default`
`<prefix>->precon->sub precon`<sup>8</sup>|string| Preconditioner used in each block, same options as for `<prefix>->precon`| `default`
`solver category` | string | as `<prefix>->category` | `default`
<sup>1</sup>Requires `dune-istl` 2.7<br>
<sup>2</sup>Requires external solver package<br>
<sup>3</sup>Only for `cg`, `pcg`, `fcg`, `cfcg`, `bcgs`, `minres`, `gmres` and `fgmres` solvers<br>
<sup>4</sup>Only for `pcg` and `gmres` solvers<br>
<sup>5</sup>Only for `superlu` solver<br>
<sup>6</sup>See `amdis/linearalgebra/istl/AMGPrecon.hpp` for subkeys and possible values<br>
<sup>7</sup>Only `solver` preconditioner, subkeys use the same values as above e.g. `<prefix>->precon->solver->max iteration: 500`<br>
<sup>8</sup>Only for `bjacobi` preconditioner, subkeys use the same values as above e.g. `<prefix>->precon->sub precon->iterations: 1`<br>
<sup>9</sup>Checks global parameter `solver category` if not set<br>
### Grids
key | type | values | default value
------------------------------|--------|----------------------------------|-------------
`<prefix>->macro file name` | string | Filename of the mesh file | None<sup>1</sup>
`<prefix>->structured`<sup>2</sup>| string | Type of the structured grid, one of:<br>`cube` (rectangular or cube grid),<br>`simplex` (triangular or tetrahedral grid) | Grid dependant
`<prefix>->global refinements`| int | Number of initial global refinement steps| `0`
`<prefix>->load balance` | bool | Perform initial load balance if `true`| No load balance
`<prefix>->min corner`<sup>2</sup> | dim-sized array of global coordinates | lower left corner of the domain | `(0,...,0)`
`<prefix>->max corner`<sup>2</sup> | dim-sized array of global coordinates | upper right corner of the domain | `(1,...,1)`
`<prefix>->num cells`<sup>2</sup> | dim-sized array of int | number of blocks in each coordinate direction | `(1,...,1)`
`<prefix>->overlap`<sup>3</sup>| int | Number of overlap elements on a distributed grid | `0`
`<prefix>->periodic`<sup>3</sup>|string| String of ones and zeros indicating if the grid is periodic in the xyz-direction, e.g. `010` periodic in y-direction | `0...0` (not periodic)
<sup>1</sup>If no file name is given a structured grid will be constructed<br>
<sup>2</sup>Only for structured grid generation<br>
<sup>3</sup>`YaspGrid` only<br>
### Adaptivity
#### Parameters used by the `AdaptInfo` class
key | type | values | default value
-----------------------------------|------|----------------------------------|-------------
`<prefix>->tolerance` |double| Tolerance for the (absolute or relative) error | `0.0`
`<prefix>->time tolerance` |double| Time tolerance | `0.0`
`<prefix>->time relative tolerance`|double| Relative time tolerance | `0.0`
`<prefix>->coarsen allowed` | int | `0` (coarsening not allowed),<br> `1` (coarsening allowed)| `0`
`<prefix>->refinement allowed` | int | `0` (refinement not allowed),<br> `1` (refinement allowed)| `0`
`<prefix>->sum factor` |double| Factor to combine max and integral time estimate|`1.0`
`<prefix>->max factor` |double| Factor to combine max and integral time estimate|`0.0`
`<prefix>->start time` |double| Initial time | `0.0`
`<prefix>->timestep` |double| Time step size to be used | `0.0`
`<prefix>->end time` |double| Final time | `0.0`
`<prefix>->max iteration` | int | Maximal allowed number of iterations of the adaptive procedure| No maximum
`<prefix>->max timestep iteration` | int | Maximal number of iterations for choosing a timestep|`30`
`<prefix>->max time iteration` | int | Maximal number of time iterations| `30`
`<prefix>->min timestep` |double| Minimal step size | Square root of floating-point epsilon
`<prefix>->max timestep` |double| Maximal step size | Square root of floating-point maximal value
`<prefix>->number of timesteps` | int | Number of fixed timestep iterations | None<sup>1</sup>
`<prefix>->time tolerance` |double| Tolerance for the overall time error| `1.0`
<sup>1</sup>If this is set to a nonzero value the computation is done `number of timesteps` times with a fixed time step size<br>
#### Parameters used by the `AdaptInstationary` class
key | type | values | default value
------------------------------|------|----------------------------------|-------------
`<prefix>->strategy` | int | `0` (explicit time strategy),<br>`1` (implicit time strategy),<br>`2` (simple adaptive time strategy) | `0`
`<prefix>->time delta 1` |double| Parameter Δ<sub>1</sub> used in time step reduction| `1/√2`
`<prefix>->time delta 2` |double| Parameter Δ<sub>2</sub> used in time step enlargement| `√2`
`<prefix>->break when stable` | bool | Stops the iteration when the instationary problem is stable when set to `true`| `false`
### Filewriters
key | type | values | default value
------------------------------|--------|----------------------------------|-------------
`<prefix>->format` | vector&lt;string&gt; | List of format strings from the following list<sup>1</sup>:<br>`vtk` (VTK format using the dune-grid writer),<br>`dune-vtk` (VTK format using the dune-vtk writer)<sup>2</sup>,<br>`gmsh` (GMsh file format),<br>`backup` (backup writer) | Not optional
`<prefix>->filename` | string | Filename of the output file excluding the file extension| `solution`
`<prefix>->output directory` | string | Directory where to put the files | `.`
`<prefix>->name` | string | Name of the data vector in the output file| `solution`
`<prefix>->write every i-th timestep`|int| Timestep number interval | `0`
`<prefix>->write after timestep`|double| Time interval | `0.0`
`<prefix>->animation`<sup>345</sup>|bool| write .pvd file if `true` | `false`
`<prefix>->mode`<sup>3</sup> | int | `0` (output in ASCII),<br>`1` (output appended base64 binary) | `0`
`<prefix>->subsampling`<sup>3</sup>|int| Number of refinement intervals for subsampling| No subsampling
`<prefix>->mode`<sup>4</sup> | int | `0` (ASCII),<br>`1` (binary),<br>`2` (compressed) | `0`
`<prefix>->precision`<sup>45</sup>| int| `0` (float),<br>`1` (double) | `0`
`<prefix>->animation`<sup>6</sup>|bool | Append timestepnumber if `true`, otherwise always override | `false`
<sup>1</sup>One filewriter will be set up for each format string in the list<br>
<sup>2</sup>Requires the `dune-vtk` module<br>
<sup>3</sup>Only for `vtk` format<br>
<sup>4</sup>Only for `dune-vtk` format<br>
<sup>5</sup>Only for `gmsh` format<br>
<sup>6</sup>Only for `backup` format<br>
......@@ -19,6 +19,7 @@ nav:
- DOFVector and DiscreteFunction: reference/DOFVector.md
- GlobalBasis: reference/GlobalBasis.md
- GridFunctions: reference/GridFunctions.md
- Initfile: reference/Initfile.md
- Operators: reference/Operators.md
- Problem: reference/Problem.md
- MatrixBase and VectorBase: reference/MatVecBase.md
......
Markdown is supported
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment