Multi-implicit SDC sweeper type, extends abstract sweeper
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | npieces | ||||
logical, | public | :: | use_LUq | ||||
real(kind=pfdp), | public, | allocatable | :: | QdiffE(:,:) | |||
real(kind=pfdp), | public, | allocatable | :: | QdiffI(:,:) | |||
real(kind=pfdp), | public, | allocatable | :: | QtilE(:,:) | |||
real(kind=pfdp), | public, | allocatable | :: | QtilI(:,:) | |||
real(kind=pfdp), | public, | allocatable | :: | dtsdc(:) | |||
class(pf_encap_t), | public, | allocatable | :: | I3(:) | |||
class(pf_encap_t), | public, | allocatable | :: | rhs |
This is the interface for the routine to compute the RHS function values Evaluate f_piece(y), where piece is one or two Evaluate f_piece(y), where piece is one or two
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_encap_t), | intent(in) | :: | y | Argument for evaluation |
||
real(kind=pfdp), | intent(in) | :: | t | Time at evaluation |
||
integer, | intent(in) | :: | level_index | Level index |
||
class(pf_encap_t), | intent(inout) | :: | f | RHS function value |
||
integer, | intent(in) | :: | piece | Which piece to evaluate |
Solve the equation where n is given by the argument piece
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_encap_t), | intent(inout) | :: | y | Solution of implicit solve |
||
real(kind=pfdp), | intent(in) | :: | t | Time of solve |
||
real(kind=pfdp), | intent(in) | :: | dtq | dt*quadrature weight |
||
class(pf_encap_t), | intent(in) | :: | rhs | RHS for solve |
||
integer, | intent(in) | :: | level_index | Level index |
||
class(pf_encap_t), | intent(inout) | :: | f | f_n of solution y |
||
integer, | intent(in) | :: | piece | Which piece to evaluate |
Assign level pointer
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
type(pf_pfasst_t), | intent(inout), | target | :: | pf | PFASST structure |
|
integer, | intent(in) | :: | level_index | which level to sweep on |
||
real(kind=pfdp), | intent(in) | :: | t0 | Time at beginning of time step |
||
real(kind=pfdp), | intent(in) | :: | dt | time step size |
||
integer, | intent(in) | :: | nsweeps | number of sweeps to do |
||
integer, | intent(in), | optional | :: | flags |
Array of substep sizes
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(inout) | :: | lev |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(in) | :: | lev | |||
class(pf_encap_t), | intent(in) | :: | qSDC(:) | |||
class(pf_encap_t), | intent(in) | :: | fSDC(:,:) | |||
real(kind=pfdp), | intent(in) | :: | dt | |||
class(pf_encap_t), | intent(inout) | :: | fintSDC(:) | |||
integer, | intent(in), | optional | :: | flags |
Subroutine to compute Residual
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(inout) | :: | lev | Current level |
||
real(kind=pfdp), | intent(in) | :: | dt | Time step |
||
integer, | intent(in), | optional | :: | flags |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(inout) | :: | lev | |||
real(kind=pfdp), | intent(in) | :: | t0 | |||
integer, | intent(in), | optional | :: | flags | ||
integer, | intent(in), | optional | :: | step |
Subroutine to evaluate the function values at all nodes
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(inout) | :: | lev | Current level |
||
real(kind=pfdp), | intent(in) | :: | t(:) | Array of times at each node |
||
integer, | intent(in), | optional | :: | flags | ||
integer, | intent(in), | optional | :: | step |
Subroutine to evaluate function value at node m
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(inout) | :: | lev | Current level |
||
real(kind=pfdp), | intent(in) | :: | t | Time at which to evaluate |
||
integer, | intent(in) | :: | m | Node at which to evaluate |
||
integer, | intent(in), | optional | :: | flags | ||
integer, | intent(in), | optional | :: | step |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(inout) | :: | lev |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_misdcQ_t), | intent(inout) | :: | this | |||
class(pf_level_t), | intent(inout) | :: | lev |
type, extends(pf_sweeper_t), abstract :: pf_misdcQ_t
real(pfdp), allocatable :: QdiffE(:,:)
real(pfdp), allocatable :: QdiffI(:,:)
real(pfdp), allocatable :: QtilE(:,:)
real(pfdp), allocatable :: QtilI(:,:)
real(pfdp), allocatable :: dtsdc(:)
class(pf_encap_t), allocatable :: I3(:)
class(pf_encap_t), allocatable :: rhs
contains
procedure(pf_f_eval_p), deferred :: f_eval
procedure(pf_f_comp_p), deferred :: f_comp
procedure :: sweep => misdcQ_sweep
procedure :: initialize => misdcQ_initialize
procedure :: integrate => misdcQ_integrate
procedure :: residual => misdcQ_residual
procedure :: spreadq0 => misdcQ_spreadq0
procedure :: evaluate_all => misdcQ_evaluate_all
procedure :: evaluate => misdcQ_evaluate
procedure :: destroy => misdcQ_destroy
procedure :: misdcQ_destroy
end type pf_misdcQ_t