Module to restrict solutions between pfasst levels and create the FAS tau correction
Restrict (in time and space) fine level to coarse and set coarse level FAS correction.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(pf_pfasst_t), | intent(inout), | target | :: | pf | ||
real(kind=pfdp), | intent(in) | :: | t0 | time at beginning of step |
||
real(kind=pfdp), | intent(in) | :: | dt | time step |
||
integer, | intent(in) | :: | level_index | defines which level to restrict |
||
integer, | intent(in), | optional | :: | flags |
Restrict (in time and space) f_sol_array to c_sol_array Depending on the flag INTEGRAL, we may be restricting solutions, or integrals of F
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_level_t), | intent(inout) | :: | f_lev_ptr | pointer to fine level |
||
class(pf_level_t), | intent(inout) | :: | c_lev_ptr | pointer to coarse level |
||
class(pf_encap_t), | intent(inout) | :: | f_encap_array(:) | array of fine level data to be restricted |
||
class(pf_encap_t), | intent(inout) | :: | c_encap_array(:) | array of coarse level data to be computed |
||
logical, | intent(in) | :: | IS_INTEGRAL | flag determines if it is integral data being restricted |
||
real(kind=pfdp), | intent(in) | :: | f_time(:) | time at the fine nodes |
||
integer, | intent(in), | optional | :: | flags |
Apply a matrix (tmat or rmat) to src and add to dst. Mathematically this is dst= dst + amatsrc Where dst and src are vectors, mat is a matrix, and a is a scalar If the optional variable "zero" is provided and is true, then we compute dst= amatsrc
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_encap_t), | intent(inout) | :: | dst(:) | destination vector |
||
real(kind=pfdp), | intent(in) | :: | a | scalar |
||
real(kind=pfdp), | intent(in) | :: | mat(:,:) | matrix |
||
class(pf_encap_t), | intent(in) | :: | src(:) | src vector |
||
logical, | intent(in), | optional | :: | zero | If false, don't zero out the the dst variable before computing |
|
integer, | intent(in), | optional | :: | flags | Used for choosing which variable to operate on |
Apply a matrix (tmat or rmat) to src and add to dst.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_encap_t), | intent(inout) | :: | dst(:) | destination vector |
||
real(kind=pfdp), | intent(in) | :: | a | scalar |
||
real(kind=pfdp), | intent(in) | :: | mat(:,:) | matrix |
||
class(pf_encap_t), | intent(in) | :: | src(:) | src vector |
||
logical, | intent(in), | optional | :: | zero | If false, don't zero out the the dst variable before computing |
|
integer, | intent(in), | optional | :: | flags | Used for choosing which variable to operate on |