|
| 1 | +!< Define the abstract type *integrand* for building FOODIE ODE integrators. |
| 2 | + |
| 3 | +module foodie_integrand_object |
| 4 | +!< Define the abstract type *integrand* for building FOODIE ODE integrators. |
| 5 | + |
| 6 | +use penf, only : I_P, R_P |
| 7 | + |
| 8 | +implicit none |
| 9 | +private |
| 10 | +public :: integrand_object |
| 11 | + |
| 12 | +type, abstract :: integrand_object |
| 13 | + !< Abstract type for building FOODIE ODE integrators. |
| 14 | +#ifdef CAF |
| 15 | + class(*), allocatable :: dummy_to_allow_extensions[:] !< Dummy member to allow concrete extensions with coarray members. |
| 16 | +#endif |
| 17 | + contains |
| 18 | + ! public deferred procedures that concrete integrand-field must implement |
| 19 | + procedure(time_derivative), pass(self), deferred, public :: t !< Time derivative, residuals. |
| 20 | + ! operators |
| 21 | + procedure(local_error_operator), pass(lhs), deferred, public :: local_error !< `||integrand - integrand||` operator. |
| 22 | + generic, public :: operator(.lterror.) => local_error !< Estimate local truncation error. |
| 23 | + ! + |
| 24 | + procedure(symmetric_operator), pass(lhs), deferred, public :: integrand_add_integrand !< `+` operator. |
| 25 | + procedure(integrand_op_real), pass(lhs), deferred, public :: integrand_add_real !< `+ real` operator. |
| 26 | + procedure(real_op_integrand), pass(rhs), deferred, public :: real_add_integrand !< `real +` operator. |
| 27 | + generic, public :: operator(+) => integrand_add_integrand, & |
| 28 | + integrand_add_real, & |
| 29 | + real_add_integrand !< Overloading `+` operator. |
| 30 | + ! * |
| 31 | + procedure(symmetric_operator), pass(lhs), deferred, public :: integrand_multiply_integrand !< `*` operator. |
| 32 | + procedure(integrand_op_real), pass(lhs), deferred, public :: integrand_multiply_real !< `* real` operator. |
| 33 | + procedure(real_op_integrand), pass(rhs), deferred, public :: real_multiply_integrand !< `real *` operator. |
| 34 | + procedure(integrand_op_real_scalar), pass(lhs), deferred, public :: integrand_multiply_real_scalar !< `* real_scalar` operator. |
| 35 | + procedure(real_scalar_op_integrand), pass(rhs), deferred, public :: real_scalar_multiply_integrand !< `real_scalar *` operator. |
| 36 | + generic, public :: operator(*) => integrand_multiply_integrand, & |
| 37 | + integrand_multiply_real, & |
| 38 | + real_multiply_integrand, & |
| 39 | + integrand_multiply_real_scalar, & |
| 40 | + real_scalar_multiply_integrand !< Overloading `*` operator. |
| 41 | + ! - |
| 42 | + procedure(symmetric_operator), pass(lhs), deferred, public :: integrand_sub_integrand !< `-` operator. |
| 43 | + procedure(integrand_op_real), pass(lhs), deferred, public :: integrand_sub_real !< `- real` operator. |
| 44 | + procedure(real_op_integrand), pass(rhs), deferred, public :: real_sub_integrand !< `real -` operator. |
| 45 | + generic, public :: operator(-) => integrand_sub_integrand, & |
| 46 | + integrand_sub_real, & |
| 47 | + real_sub_integrand !< Overloading `-` operator. |
| 48 | + ! = |
| 49 | + procedure(assignment_integrand), pass(lhs), deferred, public :: assign_integrand !< `=` operator. |
| 50 | + procedure(assignment_real), pass(lhs), deferred, public :: assign_real !< `= real` operator. |
| 51 | + generic, public :: assignment(=) => assign_integrand, assign_real !< Overloading `=` assignament. |
| 52 | +endtype integrand_object |
| 53 | + |
| 54 | +abstract interface |
| 55 | + !< Abstract type bound procedures necessary for implementing a concrete extension of [[integrand_object]]. |
| 56 | + |
| 57 | + function time_derivative(self, t) result(dState_dt) |
| 58 | + !< Time derivative function of integrand class, i.e. the residuals function. |
| 59 | + import :: integrand_object, R_P |
| 60 | + class(integrand_object), intent(in) :: self !< Integrand field. |
| 61 | + real(R_P), optional, intent(in) :: t !< Time. |
| 62 | + real(R_P), allocatable :: dState_dt(:) !< Result of the time derivative function of integrand field. |
| 63 | + endfunction time_derivative |
| 64 | + |
| 65 | + ! operators |
| 66 | + function local_error_operator(lhs, rhs) result(error) |
| 67 | + !< Estimate local truncation error between 2 solution approximations. |
| 68 | + import :: integrand_object, R_P |
| 69 | + class(integrand_object), intent(in) :: lhs !< Left hand side. |
| 70 | + class(integrand_object), intent(in) :: rhs !< Right hand side. |
| 71 | + real(R_P) :: error !< Error estimation. |
| 72 | + endfunction local_error_operator |
| 73 | + |
| 74 | + pure function integrand_op_real(lhs, rhs) result(operator_result) |
| 75 | + !< Asymmetric type operator `integrand.op.real`. |
| 76 | + import :: integrand_object, R_P |
| 77 | + class(integrand_object), intent(in) :: lhs !< Left hand side. |
| 78 | + real(R_P), intent(in) :: rhs(1:) !< Right hand side. |
| 79 | + real(R_P), allocatable :: operator_result(:) !< Operator result. |
| 80 | + endfunction integrand_op_real |
| 81 | + |
| 82 | + pure function real_op_integrand(lhs, rhs) result(operator_result) |
| 83 | + !< Asymmetric type operator `real.op.integrand`. |
| 84 | + import :: integrand_object, R_P |
| 85 | + class(integrand_object), intent(in) :: rhs !< Right hand side. |
| 86 | + real(R_P), intent(in) :: lhs(1:) !< Left hand side. |
| 87 | + real(R_P), allocatable :: operator_result(:) !< Operator result. |
| 88 | + endfunction real_op_integrand |
| 89 | + |
| 90 | + pure function integrand_op_real_scalar(lhs, rhs) result(operator_result) |
| 91 | + !< Asymmetric type operator `integrand.op.real`. |
| 92 | + import :: integrand_object, R_P |
| 93 | + class(integrand_object), intent(in) :: lhs !< Left hand side. |
| 94 | + real(R_P), intent(in) :: rhs !< Right hand side. |
| 95 | + real(R_P), allocatable :: operator_result(:) !< Operator result. |
| 96 | + endfunction integrand_op_real_scalar |
| 97 | + |
| 98 | + pure function real_scalar_op_integrand(lhs, rhs) result(operator_result) |
| 99 | + !< Asymmetric type operator `real.op.integrand`. |
| 100 | + import :: integrand_object, R_P |
| 101 | + real(R_P), intent(in) :: lhs !< Left hand side. |
| 102 | + class(integrand_object), intent(in) :: rhs !< Right hand side. |
| 103 | + real(R_P), allocatable :: operator_result(:) !< Operator result. |
| 104 | + endfunction real_scalar_op_integrand |
| 105 | + |
| 106 | + pure function symmetric_operator(lhs, rhs) result(operator_result) |
| 107 | + !< Symmetric type operator integrand.op.integrand. |
| 108 | + import :: integrand_object, R_P |
| 109 | + class(integrand_object), intent(in) :: lhs !< Left hand side. |
| 110 | + class(integrand_object), intent(in) :: rhs !< Right hand side. |
| 111 | + real(R_P), allocatable :: operator_result(:) !< Operator result. |
| 112 | + endfunction symmetric_operator |
| 113 | + |
| 114 | + pure subroutine assignment_integrand(lhs, rhs) |
| 115 | + !< Symmetric assignment integrand = integrand. |
| 116 | + import :: integrand_object |
| 117 | + class(integrand_object), intent(inout) :: lhs !< Left hand side. |
| 118 | + class(integrand_object), intent(in) :: rhs !< Right hand side. |
| 119 | + endsubroutine assignment_integrand |
| 120 | + |
| 121 | + pure subroutine assignment_real(lhs, rhs) |
| 122 | + !< Symmetric assignment integrand = integrand. |
| 123 | + import :: integrand_object, R_P |
| 124 | + class(integrand_object), intent(inout) :: lhs !< Left hand side. |
| 125 | + real(R_P), intent(in) :: rhs(1:) !< Right hand side. |
| 126 | + endsubroutine assignment_real |
| 127 | +endinterface |
| 128 | +endmodule foodie_integrand_object |
0 commit comments