Time-Dependent Solvers

TimeControlSolver

The structure TimeControlSolver can be used to setup a time-dependent solver that can be configured in a similar manner as the time-independent ones (subiterations, nonlinear iterations, linear solvers). The following table lists the available TimeIntegrationRules:

Time integration ruleorderFormula
BackwardEuler1$(M^{n+1} + A^{n+1}) u^{n+1} = F^{n+1} + M^{n+1} u^n$
CrankNicolson2$(M^{n+1} + 0.5 A^{n+1}) u^{n+1} = 0.5 (F^{n+1} + F^n - A^n u^n) + M^{n+1} u^n$

Note that currently the time-derivative (M terms) is added by the TimeControlSolver in each integration step and is in general not part of the PDEDescription (this might change in future). The default time derivative is a scaled (depends on the integration rule) mass matrix of the used finite element space, but the user can overwrite it via optional constructor arguments (experimental).

In case of the Crank-Nicolson scheme the user can mask unknowns of the PDE as an algebraic constraint (see add_unknown! in ProblemDescription). For these variables old iterates are not used on the right-hand side of the iteration formula. The pressure in the Navier-Stokes system is an example for such a constraint.

GradientRobustMultiPhysics.TimeControlSolverType
function TimeControlSolver(
    PDE::PDEDescription,                                        # PDE system description
    InitialValues::FEVector{T,Tv,Ti},                           # contains initial values and stores solution of advance! methods
    TIR::Type{<:AbstractTimeIntegrationRule} = BackwardEuler;   # Time integration rule
    dt_operator = [],                                           # Operator in time derivative (default: Identity) for each subiteration (applied to test and ansatz function)
    dt_action = [],                                             # Action in time derivative (deulta: NoAction) for each subiteration
    dt_is_nonlinear = [],                                          # is the time derivative nonlinear?
    T_time = Float64,                                           # Type for timestep and total time
    kwargs...) where {T,Tv,Ti}                                  # additional solver arguments

Creates a time-dependent solver that can be advanced in time with advance!. The FEVector Solution stores the initial state but also the solution at the current time. The argument TIR carries the time integration rule to be use (e.g. BackwardEuler or CrankNicolson). T_time determines the NumberType of the timesteps and total time.

Keyword arguments:

  • subiterations: an array of equation subsets (each an array) that should be solved together in each fixpoint iteration. Default: ''auto''

  • showiterationdetails: show details (residuals etc.) of each iteration. Default: true

  • timedependent_equations: array of the equations that should get a time derivative (only for TimeControlSolver). Default: Any[]

  • reltol: reltol for linear solver (if iterative). Default: 1.0e-12

  • show_statistics: show some statistics like assembly times. Default: false

  • skipupdate: matrix update (for the j-th sub-iteration) will be performed each skipupdate[j] iteration; -1 means only in the first iteration. Default: [1]

  • linsolver: any AbstractFactorization from LinearSolve.jl (default = UMFPACKFactorization). Default: LinearSolve.UMFPACKFactorization(true, true)

  • time: time at which time-dependent data functions are evaluated or initial time for TimeControlSolver. Default: 0

  • parallel_storage: assemble storaged operators in parallel for loop. Default: false

  • showsolverconfig: show the complete solver configuration before starting to solve. Default: false

  • abstol: abstol for linear solver (if iterative). Default: 1.0e-12

  • checknonlinearresidual: check the nonlinear residual in last nonlinear iteration (causes one more reassembly of nonlinear terms). Default: ''auto''

  • fixed_penalty: penalty that is used for the values of fixed degrees of freedom (e.g. by Dirichlet boundary data or global constraints). Default: 1.0e60

  • target_residual: stop fixpoint iterations if the absolute (nonlinear) residual is smaller than this number. Default: 1.0e-12

  • maxiterations: maximal number of nonlinear iterations (TimeControlSolver runs that many in each time step). Default: ''auto''

Further (very experimental) optional arguments for TimeControlSolver are:

  • dt_operator : (array of) operators applied to testfunctions in time derivative (default: Identity)
  • dt_action : (array of) actions that are applied to the ansatz function in the time derivative (to include parameters etc.)
  • dtisnonlinear : (array of) booleans to decide which time derivatives should be recomputed in each timestep
source
GradientRobustMultiPhysics.advance!Function
function TimeControlSolver(
    advance!(TCS::TimeControlSolver, timestep::Real = 1e-1)

Advances a TimeControlSolver one step in time with the given timestep.

source

Advancing a TimeControlSolver

There are two functions that advance the TimeControlSolver automatically until a given final time (advance_until_time!) is reached or until stationarity is reached (advance_until_stationarity!). As an experimental feature, one can add the module DifferentialEquations.jl as the first argument to these methods to let this module run the time integration (the native TimeIntegrationRule argument in the TimeControlSolver constuctor is ignored in this case).

GradientRobustMultiPhysics.advance_until_time!Function
advance_until_time!(TCS::TimeControlSolver, timestep, finaltime; finaltime_tolerance = 1e-15, do_after_each_timestep = nothing)

Advances a TimeControlSolver in time with the given (initial) timestep until the specified finaltime is reached (up to the specified tolerance). The function doaftertimestep is called after each timestep and can be used to print/save data (and maybe timestep control in future).

source
advance_until_time!(DiffEQ::Module, TCS::TimeControlSolver, timestep, finaltime; solver = nothing, abstol = 1e-1, reltol = 1e-1, dtmin = 0, adaptive::Bool = true)

Advances a TimeControlSolver in time with the given (initial) timestep until the specified finaltime is reached (up to the specified tolerance) with the given exterior time integration module. The only valid Module here is DifferentialEquations.jl and the optional arguments are passed to it. If solver == nothing the solver Rosenbrock23(autodiff = false) will be chosen. For more choices please consult the documentation of DifferentialEquations.jl.

Also note that this is a highly experimental feature and will not work for general TimeControlSolvers configuration (e.g. in the case of several subiterations or, it seems, saddle point problems). Also have a look at corressponding the example in the advanced examples section.

source
GradientRobustMultiPhysics.advance_until_stationarity!Function
advance_until_stationarity!(TCS::TimeControlSolver, timestep; stationarity_threshold = 1e-11, maxtimesteps = 100, do_after_each_timestep = nothing)

Advances a TimeControlSolver in time with the given (initial) timestep until stationarity is detected (change of variables below threshold) or a maximal number of time steps is exceeded. The function doaftertimestep is called after each timestep and can be used to print/save data (and maybe timestep control in future).

source