mform/m_mform_linalg.f90
| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | !> The m-form linear algebra module | ||
| 2 | !> | ||
| 3 | !> This module provides an inner product, L2 projection, and L2-error computation for the m-form spaces | ||
| 4 | module m_mform_linalg | ||
| 5 | use m_common, only: wp, user_function_3d_interface | ||
| 6 | use m_bspline_basis, only: BSplineFun | ||
| 7 | use m_mform_basis, only: MFormFun, MFormSpace | ||
| 8 | #include "petsc.fi" | ||
| 9 | |||
| 10 | implicit none | ||
| 11 | |||
| 12 | private | ||
| 13 | |||
| 14 | public :: inner_product, l2_error, l2_norm, l2_projection | ||
| 15 | public :: solve_curlA_equals_B | ||
| 16 | |||
| 17 | !> @brief Compute the (L2) inner product of all elements of a m-form space with a user-defined function | ||
| 18 | interface inner_product | ||
| 19 | procedure inner_product_mform_scalar, inner_product_mform_vector | ||
| 20 | end interface | ||
| 21 | |||
| 22 | interface l2_projection | ||
| 23 | procedure l2_projection_mform_scalar, l2_projection_mform_vector | ||
| 24 | end interface | ||
| 25 | |||
| 26 | !> @brief Compute the L2 norm of the difference between an m-form and a user-defined function | ||
| 27 | interface l2_error | ||
| 28 | procedure l2_error_mform_scalar, l2_error_mform_vector | ||
| 29 | end interface | ||
| 30 | |||
| 31 | !> @brief Compute the L2 norm of an m-form | ||
| 32 | interface l2_norm | ||
| 33 | procedure l2_norm_mform | ||
| 34 | end interface | ||
| 35 | |||
| 36 | interface solve_curlA_equals_B | ||
| 37 | module procedure solve_curlA_equals_B_tensorprod | ||
| 38 | end interface | ||
| 39 | |||
| 40 | ! NOTE: a compiler bug in NVHPCSDK requires these variables to be in the module scope | ||
| 41 | ! NOTE: see also: https://forums.developer.nvidia.com/t/compiler-bug-passing-a-procedure-with-a-variable-from-a-local-scope-as-input-to-a-function-leads-to-unexpected-behavior/359876 | ||
| 42 | type(BSplineFun) :: twoform_x_integrated, twoform_y_copy | ||
| 43 | |||
| 44 | contains | ||
| 45 | |||
| 46 | !> @brief Compute the (L2) inner product of all elements of a scalar-valued m-form with a user-defined function | ||
| 47 | !> | ||
| 48 | !> @param[inout] v The resulting m-form containing inner products values | ||
| 49 | !> @param[in] space The m-form space (defines the first argument of the inner product) | ||
| 50 | !> @param[in] userfun The second argument of the inner product is a user-defined function | ||
| 51 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 52 | !> @param[in] n_quad_extra _(optional)_ The extra number of quadrature points to use for the inner product on top of those needed | ||
| 53 | !> to exactly integrate the B-spline spaces | ||
| 54 | !> @param[in] user_fun_is_physical _(optional)_ Whether the user function is defined in physical coordinates (default: .false.) | ||
| 55 | !> | ||
| 56 | !> @note The resulting m-form `v` is initialized only for the response B-spline indices of the tensor product space (if needed, | ||
| 57 | !> call v%distribute() to distribute the data) | ||
| 58 | 433 | subroutine inner_product_mform_scalar(v, space, userfun, coord_transform, n_quad_extra, user_fun_is_physical) | |
| 59 | use m_tensorprod_linalg, only: inner_product | ||
| 60 | use m_tensorprod_quadrature, only: DEFAULT_EXTRA_QUADRATURE_POINTS | ||
| 61 | use m_coord_transform, only: CoordTransformAbstract, jacobian, transform | ||
| 62 | implicit none | ||
| 63 | |||
| 64 | type(MFormFun), intent(inout) :: v | ||
| 65 | type(MFormSpace), intent(in) :: space | ||
| 66 | procedure(user_function_3d_interface) :: userfun | ||
| 67 | class(CoordTransformAbstract), intent(in), optional :: coord_transform | ||
| 68 | integer, intent(in), optional :: n_quad_extra | ||
| 69 | logical, intent(in), optional :: user_fun_is_physical | ||
| 70 | |||
| 71 | logical :: user_fun_is_physical_ | ||
| 72 | integer :: n_quad_extra_ | ||
| 73 | |||
| 74 |
2/2✓ Branch 2 → 3 taken 76 times.
✓ Branch 2 → 5 taken 357 times.
|
433 | if (v%is_initialized()) then |
| 75 | ! NOTE: when the m-form function is already initialized, it is ASSUMED that its space matches with the input space | ||
| 76 | ! TODO: enforce this? we need a cheap yet robust way to check this | ||
| 77 | 76 | call v%reset() | |
| 78 | else | ||
| 79 | 357 | call v%init(space) | |
| 80 | end if | ||
| 81 | |||
| 82 | 433 | user_fun_is_physical_ = .false. | |
| 83 |
5/6✓ Branch 7 → 8 taken 77 times.
✓ Branch 7 → 11 taken 356 times.
✓ Branch 8 → 9 taken 77 times.
✗ Branch 8 → 11 not taken.
✓ Branch 9 → 10 taken 4 times.
✓ Branch 9 → 11 taken 73 times.
|
433 | if (present(coord_transform) .and. present(user_fun_is_physical)) then |
| 84 | 4 | user_fun_is_physical_ = user_fun_is_physical | |
| 85 | end if | ||
| 86 | |||
| 87 |
1/2✗ Branch 11 → 12 not taken.
✓ Branch 11 → 18 taken 433 times.
|
433 | if (space%m == 1 .or. space%m == 2) then |
| 88 | write (*, '(A,I0,A)') "ERROR in inner_product_mform_scalar: this subroutine is for scalar-valued forms. " & | ||
| 89 | ✗ | //"The provided space is a ", space%m, "-form, which is vector-valued and requires 3 user functions. " | |
| 90 | ✗ | error stop 1 | |
| 91 | end if | ||
| 92 | |||
| 93 |
2/2✓ Branch 18 → 19 taken 38 times.
✓ Branch 18 → 20 taken 395 times.
|
433 | if (present(n_quad_extra)) then |
| 94 | 38 | n_quad_extra_ = n_quad_extra | |
| 95 | else | ||
| 96 |
3/4✓ Branch 20 → 21 taken 39 times.
✓ Branch 20 → 23 taken 356 times.
✓ Branch 21 → 22 taken 39 times.
✗ Branch 21 → 23 not taken.
|
395 | if (present(coord_transform)) then |
| 97 | 39 | n_quad_extra_ = DEFAULT_EXTRA_QUADRATURE_POINTS | |
| 98 | else | ||
| 99 | 356 | n_quad_extra_ = 0 | |
| 100 | end if | ||
| 101 | end if | ||
| 102 | |||
| 103 | 433 | call inner_product(v%tp_funs(1), space%tp_spaces(1), wrapped_userfun, n_quad_extra=n_quad_extra_) | |
| 104 | |||
| 105 | contains | ||
| 106 | 24184104 | pure real(wp) function wrapped_userfun(xp, yp, zp) result(ans) | |
| 107 | real(wp), intent(in) :: xp, yp, zp | ||
| 108 | |||
| 109 | real(wp) :: x, y, z | ||
| 110 | |||
| 111 |
3/4✓ Branch 2 → 3 taken 19606920 times.
✓ Branch 2 → 4 taken 4577184 times.
✗ Branch 3 → 4 not taken.
✓ Branch 3 → 5 taken 19606920 times.
|
24184104 | if (.not. present(coord_transform)) then |
| 112 | 4577184 | ans = userfun(xp, yp, zp) | |
| 113 | else | ||
| 114 |
2/2✓ Branch 5 → 6 taken 506880 times.
✓ Branch 5 → 7 taken 19100040 times.
|
19606920 | if (user_fun_is_physical_) then |
| 115 | 506880 | x = transform(coord_transform, 1, xp, yp, zp) | |
| 116 | 506880 | y = transform(coord_transform, 2, xp, yp, zp) | |
| 117 | 506880 | z = transform(coord_transform, 3, xp, yp, zp) | |
| 118 | 506880 | ans = userfun(x, y, z) | |
| 119 | else | ||
| 120 | 19100040 | ans = userfun(xp, yp, zp) | |
| 121 |
2/2✓ Branch 8 → 9 taken 6942900 times.
✓ Branch 8 → 10 taken 12157140 times.
|
19100040 | if (v%space%m == 3) then |
| 122 | ! The pushforward transformation for 3-forms yields a factor 1 / jacobian | ||
| 123 | 6942900 | ans = ans / jacobian(coord_transform, xp, yp, zp) | |
| 124 | end if | ||
| 125 | end if | ||
| 126 |
2/2✓ Branch 10 → 11 taken 12479700 times.
✓ Branch 10 → 12 taken 7127220 times.
|
19606920 | if (v%space%m == 0) then |
| 127 | 12479700 | ans = ans * jacobian(coord_transform, xp, yp, zp) | |
| 128 | end if | ||
| 129 | end if | ||
| 130 | |||
| 131 | 24184104 | end function wrapped_userfun | |
| 132 | end subroutine | ||
| 133 | |||
| 134 | !> @brief Compute the (L2) inner product of all elements of a vector-valued m-form with user-defined functions | ||
| 135 | !> | ||
| 136 | !> @param[inout] v The resulting m-form containing inner products values | ||
| 137 | !> @param[in] space The m-form space (defines the first argument of the inner product) | ||
| 138 | !> @param[in] userfun_x The x-component of the second argument of the inner product | ||
| 139 | !> @param[in] userfun_y The y-component of the second argument of the inner product | ||
| 140 | !> @param[in] userfun_z The z-component of the second argument of the inner product | ||
| 141 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 142 | !> @param[in] n_quad_extra _(optional)_ The extra number of quadrature points to use for the inner product on top of those needed | ||
| 143 | !> to exactly integrate the B-spline spaces | ||
| 144 | !> @param[in] user_fun_is_physical _(optional)_ Whether the user functions are defined in physical coordinates (default: .false.) | ||
| 145 | !> | ||
| 146 | !> @note The resulting m-form `v` is initialized only for the response B-spline indices of the tensor product space (if needed, | ||
| 147 | !> call v%distribute() to distribute the data) | ||
| 148 | 800 | subroutine inner_product_mform_vector(v, space, userfun_x, userfun_y, userfun_z, coord_transform, n_quad_extra, & | |
| 149 | user_fun_is_physical) | ||
| 150 | use m_tensorprod_linalg, only: inner_product | ||
| 151 | use m_tensorprod_quadrature, only: DEFAULT_EXTRA_QUADRATURE_POINTS | ||
| 152 | use m_coord_transform, only: CoordTransformAbstract, jacobian, jacobian_matrix, jacobian_matrix_inv, transform, G_matrix, & | ||
| 153 | G_matrix_inv | ||
| 154 | implicit none | ||
| 155 | |||
| 156 | type(MFormFun), intent(inout) :: v | ||
| 157 | type(MFormSpace), intent(in) :: space | ||
| 158 | procedure(user_function_3d_interface) :: userfun_x, userfun_y, userfun_z | ||
| 159 | class(CoordTransformAbstract), intent(in), optional :: coord_transform | ||
| 160 | integer, intent(in), optional :: n_quad_extra | ||
| 161 | logical, intent(in), optional :: user_fun_is_physical | ||
| 162 | |||
| 163 | logical :: user_fun_is_physical_, is_diagonal | ||
| 164 | integer :: n_quad_extra_ | ||
| 165 | |||
| 166 |
2/2✓ Branch 2 → 3 taken 49 times.
✓ Branch 2 → 5 taken 751 times.
|
800 | if (v%is_initialized()) then |
| 167 | 49 | call v%reset() | |
| 168 | else | ||
| 169 | 751 | call v%init(space) | |
| 170 | end if | ||
| 171 | |||
| 172 | 800 | is_diagonal = .true. | |
| 173 | user_fun_is_physical_ = .false. | ||
| 174 |
3/4✓ Branch 7 → 8 taken 542 times.
✓ Branch 7 → 12 taken 258 times.
✓ Branch 8 → 9 taken 542 times.
✗ Branch 8 → 12 not taken.
|
800 | if (present(coord_transform)) then |
| 175 |
2/2✓ Branch 9 → 10 taken 4 times.
✓ Branch 9 → 11 taken 538 times.
|
542 | if (present(user_fun_is_physical)) then |
| 176 | 4 | user_fun_is_physical_ = user_fun_is_physical | |
| 177 | end if | ||
| 178 | 542 | is_diagonal = coord_transform%is_orthogonal | |
| 179 | end if | ||
| 180 | |||
| 181 |
1/2✗ Branch 12 → 13 not taken.
✓ Branch 12 → 19 taken 800 times.
|
800 | if (space%m == 0 .or. space%m == 3) then |
| 182 | write (*, '(A,I0,A)') "ERROR in inner_product_mform_vector: this subroutine is for vector-valued forms. " & | ||
| 183 | ✗ | //"The provided space is a ", space%m, "-form, which is scalar-valued and requires only 1 user function. " | |
| 184 | ✗ | error stop 1 | |
| 185 | end if | ||
| 186 | |||
| 187 |
2/2✓ Branch 19 → 20 taken 68 times.
✓ Branch 19 → 21 taken 732 times.
|
800 | if (present(n_quad_extra)) then |
| 188 | 68 | n_quad_extra_ = n_quad_extra | |
| 189 | else | ||
| 190 |
3/4✓ Branch 21 → 22 taken 474 times.
✓ Branch 21 → 24 taken 258 times.
✓ Branch 22 → 23 taken 474 times.
✗ Branch 22 → 24 not taken.
|
732 | if (present(coord_transform)) then |
| 191 | 474 | n_quad_extra_ = DEFAULT_EXTRA_QUADRATURE_POINTS | |
| 192 | else | ||
| 193 | 258 | n_quad_extra_ = 0 | |
| 194 | end if | ||
| 195 | end if | ||
| 196 |
2/2✓ Branch 25 → 26 taken 2 times.
✓ Branch 25 → 29 taken 798 times.
|
800 | if (user_fun_is_physical_) then |
| 197 | 2 | call inner_product(v%tp_funs(1), space%tp_spaces(1), wrapped_physical_userfun_xx, n_quad_extra=n_quad_extra_) | |
| 198 | 2 | call inner_product(v%tp_funs(2), space%tp_spaces(2), wrapped_physical_userfun_yy, n_quad_extra=n_quad_extra_) | |
| 199 | 2 | call inner_product(v%tp_funs(3), space%tp_spaces(3), wrapped_physical_userfun_zz, n_quad_extra=n_quad_extra_) | |
| 200 | else | ||
| 201 | 798 | call inner_product(v%tp_funs(1), space%tp_spaces(1), wrapped_userfun_xx, n_quad_extra=n_quad_extra_) | |
| 202 | 798 | call inner_product(v%tp_funs(2), space%tp_spaces(2), wrapped_userfun_yy, n_quad_extra=n_quad_extra_) | |
| 203 | 798 | call inner_product(v%tp_funs(3), space%tp_spaces(3), wrapped_userfun_zz, n_quad_extra=n_quad_extra_) | |
| 204 | end if | ||
| 205 | contains | ||
| 206 | 64231318 | pure real(wp) function wrapped_userfun_xx(xp, yp, zp) result(ans) | |
| 207 | implicit none | ||
| 208 | |||
| 209 | real(wp), intent(in) :: xp, yp, zp | ||
| 210 | |||
| 211 | 64231318 | ans = userfun_x(xp, yp, zp) | |
| 212 |
3/4✓ Branch 3 → 4 taken 60653118 times.
✓ Branch 3 → 16 taken 3578200 times.
✓ Branch 4 → 5 taken 60653118 times.
✗ Branch 4 → 16 not taken.
|
64231318 | if (present(coord_transform)) then |
| 213 |
2/2✓ Branch 5 → 6 taken 50227590 times.
✓ Branch 5 → 11 taken 10425528 times.
|
60653118 | if (space%m == 1) then |
| 214 | ! The pushforward transformation for 1-forms yields | ||
| 215 | 50227590 | ans = ans * G_matrix_inv(coord_transform, 1, 1, xp, yp, zp) | |
| 216 |
2/2✓ Branch 6 → 7 taken 3685536 times.
✓ Branch 6 → 10 taken 46542054 times.
|
50227590 | if (.not. is_diagonal) then |
| 217 | 3685536 | ans = ans + userfun_y(xp, yp, zp) * G_matrix_inv(coord_transform, 1, 2, xp, yp, zp) | |
| 218 | 3685536 | ans = ans + userfun_z(xp, yp, zp) * G_matrix_inv(coord_transform, 1, 3, xp, yp, zp) | |
| 219 | end if | ||
| 220 | 50227590 | ans = ans * jacobian(coord_transform, xp, yp, zp) | |
| 221 | else | ||
| 222 | ! The pushforward transformation for 2-forms yields | ||
| 223 | 10425528 | ans = ans * G_matrix(coord_transform, 1, 1, xp, yp, zp) | |
| 224 |
2/2✓ Branch 11 → 12 taken 1788624 times.
✓ Branch 11 → 15 taken 8636904 times.
|
10425528 | if (.not. is_diagonal) then |
| 225 | 1788624 | ans = ans + userfun_y(xp, yp, zp) * G_matrix(coord_transform, 1, 2, xp, yp, zp) | |
| 226 | 1788624 | ans = ans + userfun_z(xp, yp, zp) * G_matrix(coord_transform, 1, 3, xp, yp, zp) | |
| 227 | end if | ||
| 228 | 10425528 | ans = ans / jacobian(coord_transform, xp, yp, zp) | |
| 229 | end if | ||
| 230 | end if | ||
| 231 | |||
| 232 | 64231318 | end function | |
| 233 | |||
| 234 | 489984 | pure real(wp) function wrapped_physical_userfun_xx(xp, yp, zp) result(ans) | |
| 235 | implicit none | ||
| 236 | |||
| 237 | real(wp), intent(in) :: xp, yp, zp | ||
| 238 | |||
| 239 | real(wp) :: x, y, z | ||
| 240 | |||
| 241 | 489984 | x = transform(coord_transform, 1, xp, yp, zp) | |
| 242 | 489984 | y = transform(coord_transform, 2, xp, yp, zp) | |
| 243 | 489984 | z = transform(coord_transform, 3, xp, yp, zp) | |
| 244 | |||
| 245 |
2/2✓ Branch 2 → 3 taken 268800 times.
✓ Branch 2 → 7 taken 221184 times.
|
489984 | if (space%m == 1) then |
| 246 | ! The pushforward transformation for 1-forms yields | ||
| 247 | 268800 | ans = userfun_x(x, y, z) * jacobian_matrix_inv(coord_transform, 1, 1, xp, yp, zp) | |
| 248 | 268800 | ans = ans + userfun_y(x, y, z) * jacobian_matrix_inv(coord_transform, 1, 2, xp, yp, zp) | |
| 249 | 268800 | ans = ans + userfun_z(x, y, z) * jacobian_matrix_inv(coord_transform, 1, 3, xp, yp, zp) | |
| 250 | 268800 | ans = ans * jacobian(coord_transform, xp, yp, zp) | |
| 251 | else | ||
| 252 | ! The pushforward transformation for 2-forms yields | ||
| 253 | 221184 | ans = userfun_x(x, y, z) * jacobian_matrix(coord_transform, 1, 1, xp, yp, zp) | |
| 254 | 221184 | ans = ans + userfun_y(x, y, z) * jacobian_matrix(coord_transform, 2, 1, xp, yp, zp) | |
| 255 | 221184 | ans = ans + userfun_z(x, y, z) * jacobian_matrix(coord_transform, 3, 1, xp, yp, zp) | |
| 256 | end if | ||
| 257 | 489984 | end function | |
| 258 | |||
| 259 | 63986868 | pure real(wp) function wrapped_userfun_yy(xp, yp, zp) result(ans) | |
| 260 | implicit none | ||
| 261 | |||
| 262 | real(wp), intent(in) :: xp, yp, zp | ||
| 263 | |||
| 264 | 63986868 | ans = userfun_y(xp, yp, zp) | |
| 265 |
3/4✓ Branch 3 → 4 taken 60393098 times.
✓ Branch 3 → 16 taken 3593770 times.
✓ Branch 4 → 5 taken 60393098 times.
✗ Branch 4 → 16 not taken.
|
63986868 | if (present(coord_transform)) then |
| 266 |
2/2✓ Branch 5 → 6 taken 49994696 times.
✓ Branch 5 → 11 taken 10398402 times.
|
60393098 | if (space%m == 1) then |
| 267 | ! The pushforward transformation for 1-forms yields | ||
| 268 | 49994696 | ans = ans * G_matrix_inv(coord_transform, 2, 2, xp, yp, zp) | |
| 269 |
2/2✓ Branch 6 → 7 taken 3685536 times.
✓ Branch 6 → 10 taken 46309160 times.
|
49994696 | if (.not. is_diagonal) then |
| 270 | 3685536 | ans = ans + userfun_x(xp, yp, zp) * G_matrix_inv(coord_transform, 2, 1, xp, yp, zp) | |
| 271 | 3685536 | ans = ans + userfun_z(xp, yp, zp) * G_matrix_inv(coord_transform, 2, 3, xp, yp, zp) | |
| 272 | end if | ||
| 273 | 49994696 | ans = ans * jacobian(coord_transform, xp, yp, zp) | |
| 274 | else | ||
| 275 | ! The pushforward transformation for 2-forms yields | ||
| 276 | 10398402 | ans = ans * G_matrix(coord_transform, 2, 2, xp, yp, zp) | |
| 277 |
2/2✓ Branch 11 → 12 taken 1788624 times.
✓ Branch 11 → 15 taken 8609778 times.
|
10398402 | if (.not. is_diagonal) then |
| 278 | 1788624 | ans = ans + userfun_x(xp, yp, zp) * G_matrix(coord_transform, 2, 1, xp, yp, zp) | |
| 279 | 1788624 | ans = ans + userfun_z(xp, yp, zp) * G_matrix(coord_transform, 2, 3, xp, yp, zp) | |
| 280 | end if | ||
| 281 | 10398402 | ans = ans / jacobian(coord_transform, xp, yp, zp) | |
| 282 | end if | ||
| 283 | end if | ||
| 284 | |||
| 285 | 63986868 | end function | |
| 286 | |||
| 287 | 491520 | pure real(wp) function wrapped_physical_userfun_yy(xp, yp, zp) result(ans) | |
| 288 | implicit none | ||
| 289 | |||
| 290 | real(wp), intent(in) :: xp, yp, zp | ||
| 291 | |||
| 292 | real(wp) :: x, y, z | ||
| 293 | |||
| 294 | 491520 | x = transform(coord_transform, 1, xp, yp, zp) | |
| 295 | 491520 | y = transform(coord_transform, 2, xp, yp, zp) | |
| 296 | 491520 | z = transform(coord_transform, 3, xp, yp, zp) | |
| 297 | |||
| 298 |
2/2✓ Branch 2 → 3 taken 276480 times.
✓ Branch 2 → 7 taken 215040 times.
|
491520 | if (space%m == 1) then |
| 299 | ! The pushforward transformation for 1-forms yields | ||
| 300 | 276480 | ans = userfun_x(x, y, z) * jacobian_matrix_inv(coord_transform, 2, 1, xp, yp, zp) | |
| 301 | 276480 | ans = ans + userfun_y(x, y, z) * jacobian_matrix_inv(coord_transform, 2, 2, xp, yp, zp) | |
| 302 | 276480 | ans = ans + userfun_z(x, y, z) * jacobian_matrix_inv(coord_transform, 2, 3, xp, yp, zp) | |
| 303 | 276480 | ans = ans * jacobian(coord_transform, xp, yp, zp) | |
| 304 | else | ||
| 305 | ! The pushforward transformation for 2-forms yields | ||
| 306 | 215040 | ans = userfun_x(x, y, z) * jacobian_matrix(coord_transform, 1, 2, xp, yp, zp) | |
| 307 | 215040 | ans = ans + userfun_y(x, y, z) * jacobian_matrix(coord_transform, 2, 2, xp, yp, zp) | |
| 308 | 215040 | ans = ans + userfun_z(x, y, z) * jacobian_matrix(coord_transform, 3, 2, xp, yp, zp) | |
| 309 | end if | ||
| 310 | |||
| 311 | 491520 | end function | |
| 312 | |||
| 313 | 43671429 | pure real(wp) function wrapped_userfun_zz(xp, yp, zp) result(ans) | |
| 314 | implicit none | ||
| 315 | |||
| 316 | real(wp), intent(in) :: xp, yp, zp | ||
| 317 | |||
| 318 | 43671429 | ans = userfun_z(xp, yp, zp) | |
| 319 |
3/4✓ Branch 3 → 4 taken 40117329 times.
✓ Branch 3 → 16 taken 3554100 times.
✓ Branch 4 → 5 taken 40117329 times.
✗ Branch 4 → 16 not taken.
|
43671429 | if (present(coord_transform)) then |
| 320 |
2/2✓ Branch 5 → 6 taken 24852588 times.
✓ Branch 5 → 11 taken 15264741 times.
|
40117329 | if (space%m == 1) then |
| 321 | ! The pushforward transformation for 1-forms yields | ||
| 322 | 24852588 | ans = ans * G_matrix_inv(coord_transform, 3, 3, xp, yp, zp) | |
| 323 |
2/2✓ Branch 6 → 7 taken 3685536 times.
✓ Branch 6 → 10 taken 21167052 times.
|
24852588 | if (.not. is_diagonal) then |
| 324 | 3685536 | ans = ans + userfun_x(xp, yp, zp) * G_matrix_inv(coord_transform, 3, 1, xp, yp, zp) | |
| 325 | 3685536 | ans = ans + userfun_y(xp, yp, zp) * G_matrix_inv(coord_transform, 3, 2, xp, yp, zp) | |
| 326 | end if | ||
| 327 | 24852588 | ans = ans * jacobian(coord_transform, xp, yp, zp) | |
| 328 | else | ||
| 329 | ! The pushforward transformation for 2-forms yields | ||
| 330 | 15264741 | ans = ans * G_matrix(coord_transform, 3, 3, xp, yp, zp) | |
| 331 |
2/2✓ Branch 11 → 12 taken 1788624 times.
✓ Branch 11 → 15 taken 13476117 times.
|
15264741 | if (.not. is_diagonal) then |
| 332 | 1788624 | ans = ans + userfun_x(xp, yp, zp) * G_matrix(coord_transform, 3, 1, xp, yp, zp) | |
| 333 | 1788624 | ans = ans + userfun_y(xp, yp, zp) * G_matrix(coord_transform, 3, 2, xp, yp, zp) | |
| 334 | end if | ||
| 335 | 15264741 | ans = ans / jacobian(coord_transform, xp, yp, zp) | |
| 336 | end if | ||
| 337 | end if | ||
| 338 | |||
| 339 | 43671429 | end function | |
| 340 | |||
| 341 | 488448 | pure real(wp) function wrapped_physical_userfun_zz(xp, yp, zp) result(ans) | |
| 342 | implicit none | ||
| 343 | |||
| 344 | real(wp), intent(in) :: xp, yp, zp | ||
| 345 | |||
| 346 | real(wp) :: x, y, z | ||
| 347 | |||
| 348 | 488448 | x = transform(coord_transform, 1, xp, yp, zp) | |
| 349 | 488448 | y = transform(coord_transform, 2, xp, yp, zp) | |
| 350 | 488448 | z = transform(coord_transform, 3, xp, yp, zp) | |
| 351 | |||
| 352 |
2/2✓ Branch 2 → 3 taken 258048 times.
✓ Branch 2 → 7 taken 230400 times.
|
488448 | if (space%m == 1) then |
| 353 | ! The pushforward transformation for 1-forms yields | ||
| 354 | 258048 | ans = userfun_x(x, y, z) * jacobian_matrix_inv(coord_transform, 3, 1, xp, yp, zp) | |
| 355 | 258048 | ans = ans + userfun_y(x, y, z) * jacobian_matrix_inv(coord_transform, 3, 2, xp, yp, zp) | |
| 356 | 258048 | ans = ans + userfun_z(x, y, z) * jacobian_matrix_inv(coord_transform, 3, 3, xp, yp, zp) | |
| 357 | 258048 | ans = ans * jacobian(coord_transform, xp, yp, zp) | |
| 358 | else | ||
| 359 | ! The pushforward transformation for 2-forms yields | ||
| 360 | 230400 | ans = userfun_x(x, y, z) * jacobian_matrix(coord_transform, 1, 3, xp, yp, zp) | |
| 361 | 230400 | ans = ans + userfun_y(x, y, z) * jacobian_matrix(coord_transform, 2, 3, xp, yp, zp) | |
| 362 | 230400 | ans = ans + userfun_z(x, y, z) * jacobian_matrix(coord_transform, 3, 3, xp, yp, zp) | |
| 363 | end if | ||
| 364 | |||
| 365 | 488448 | end function | |
| 366 | |||
| 367 | end subroutine | ||
| 368 | |||
| 369 | !> @brief Compute the L2 projection of a user-defined function onto a scalar-valued m-form space | ||
| 370 | !> | ||
| 371 | !> @param[out] c The resulting m-form containing the L2 projection | ||
| 372 | !> @param[in] space The m-form space onto which the projection is computed | ||
| 373 | !> @param[in] userfun The user-defined function to project onto the m-form space | ||
| 374 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 375 | !> @param[in] solver _(optional)_ The solver to use for the L2 projection, if not provided, a solver is created internally | ||
| 376 | !> @param[in] n_quad_extra _(optional)_ The extra number of quadrature points to use for the L2 projection on top of those needed | ||
| 377 | !> to exactly integrate the B-spline spaces | ||
| 378 | !> @param[in] user_fun_is_physical _(optional)_ Whether the user function is defined in physical coordinates (default: .false.) | ||
| 379 |
2/2✓ Branch 4 → 5 taken 244 times.
✓ Branch 4 → 6 taken 35 times.
|
279 | subroutine l2_projection_mform_scalar(c, space, userfun, coord_transform, solver, n_quad_extra, user_fun_is_physical) |
| 380 | use m_tensorprod_linalg, only: l2_projection | ||
| 381 | use m_mform_solver, only: GenericSolver | ||
| 382 | use m_mform_matrix, only: MassMatrix | ||
| 383 | use m_coord_transform_abstract, only: CoordTransformAbstract | ||
| 384 | implicit none | ||
| 385 | |||
| 386 | type(MFormFun), intent(out) :: c | ||
| 387 | type(MFormSpace), intent(in) :: space | ||
| 388 | procedure(user_function_3d_interface) :: userfun | ||
| 389 | |||
| 390 | class(CoordTransformAbstract), intent(in), optional :: coord_transform | ||
| 391 | type(GenericSolver), intent(in), optional :: solver | ||
| 392 | integer, intent(in), optional :: n_quad_extra | ||
| 393 | logical, intent(in), optional :: user_fun_is_physical | ||
| 394 | |||
| 395 | 279 | type(MFormFun) :: f | |
| 396 | 279 | type(MassMatrix) :: mass | |
| 397 |
2/2✓ Branch 2 → 3 taken 263 times.
✓ Branch 2 → 4 taken 16 times.
|
279 | type(GenericSolver) :: solver_ |
| 398 | |||
| 399 |
1/2✗ Branch 6 → 7 not taken.
✓ Branch 6 → 8 taken 279 times.
|
279 | if (space%m == 1 .or. space%m == 2) then |
| 400 | ✗ | error stop "l2_projection for vector-valued forms requires 3 userfunctions" | |
| 401 | end if | ||
| 402 | |||
| 403 | call inner_product(f, space, userfun, n_quad_extra=n_quad_extra, user_fun_is_physical=user_fun_is_physical, & | ||
| 404 | 279 | coord_transform=coord_transform) | |
| 405 | |||
| 406 |
1/2✓ Branch 9 → 10 taken 279 times.
✗ Branch 9 → 12 not taken.
|
279 | if (present(solver)) then |
| 407 | 279 | call solver%solve(c, f) | |
| 408 | else | ||
| 409 | ✗ | call mass%init(space, n_quad_extra=n_quad_extra, coord_transform=coord_transform) | |
| 410 | ✗ | call solver_%init(mass, coord_transform=coord_transform) | |
| 411 | ✗ | call solver_%solve(c, f) | |
| 412 | |||
| 413 | ✗ | call mass%destroy() | |
| 414 | ✗ | call solver_%destroy() | |
| 415 | end if | ||
| 416 | |||
| 417 | 279 | call f%destroy() | |
| 418 | |||
| 419 |
7/80✗ Branch 19 → 20 not taken.
✓ Branch 19 → 21 taken 279 times.
✗ Branch 21 → 22 not taken.
✓ Branch 21 → 23 taken 279 times.
✗ Branch 23 → 24 not taken.
✓ Branch 23 → 95 taken 279 times.
✗ Branch 24 → 25 not taken.
✗ Branch 24 → 26 not taken.
✗ Branch 27 → 28 not taken.
✗ Branch 27 → 71 not taken.
✗ Branch 28 → 29 not taken.
✗ Branch 28 → 31 not taken.
✗ Branch 29 → 30 not taken.
✗ Branch 29 → 31 not taken.
✗ Branch 31 → 32 not taken.
✗ Branch 31 → 70 not taken.
✗ Branch 32 → 33 not taken.
✗ Branch 32 → 34 not taken.
✗ Branch 34 → 35 not taken.
✗ Branch 34 → 47 not taken.
✗ Branch 36 → 37 not taken.
✗ Branch 36 → 46 not taken.
✗ Branch 37 → 38 not taken.
✗ Branch 37 → 39 not taken.
✗ Branch 39 → 40 not taken.
✗ Branch 39 → 41 not taken.
✗ Branch 41 → 42 not taken.
✗ Branch 41 → 43 not taken.
✗ Branch 43 → 44 not taken.
✗ Branch 43 → 45 not taken.
✗ Branch 47 → 48 not taken.
✗ Branch 47 → 58 not taken.
✗ Branch 49 → 50 not taken.
✗ Branch 49 → 57 not taken.
✗ Branch 50 → 51 not taken.
✗ Branch 50 → 52 not taken.
✗ Branch 52 → 53 not taken.
✗ Branch 52 → 54 not taken.
✗ Branch 54 → 55 not taken.
✗ Branch 54 → 56 not taken.
✗ Branch 58 → 59 not taken.
✗ Branch 58 → 69 not taken.
✗ Branch 60 → 61 not taken.
✗ Branch 60 → 68 not taken.
✗ Branch 61 → 62 not taken.
✗ Branch 61 → 63 not taken.
✗ Branch 63 → 64 not taken.
✗ Branch 63 → 65 not taken.
✗ Branch 65 → 66 not taken.
✗ Branch 65 → 67 not taken.
✗ Branch 71 → 72 not taken.
✗ Branch 71 → 73 not taken.
✗ Branch 73 → 74 not taken.
✗ Branch 73 → 75 not taken.
✗ Branch 75 → 76 not taken.
✗ Branch 75 → 77 not taken.
✗ Branch 77 → 78 not taken.
✗ Branch 77 → 79 not taken.
✗ Branch 79 → 80 not taken.
✗ Branch 79 → 81 not taken.
✗ Branch 81 → 82 not taken.
✗ Branch 81 → 94 not taken.
✗ Branch 83 → 84 not taken.
✗ Branch 83 → 93 not taken.
✗ Branch 84 → 85 not taken.
✗ Branch 84 → 86 not taken.
✗ Branch 86 → 87 not taken.
✗ Branch 86 → 88 not taken.
✗ Branch 88 → 89 not taken.
✗ Branch 88 → 90 not taken.
✗ Branch 90 → 91 not taken.
✗ Branch 90 → 92 not taken.
✗ Branch 95 → 96 not taken.
✓ Branch 95 → 97 taken 279 times.
✗ Branch 97 → 98 not taken.
✓ Branch 97 → 99 taken 279 times.
✓ Branch 99 → 100 taken 279 times.
✗ Branch 99 → 101 not taken.
✗ Branch 101 → 102 not taken.
✓ Branch 101 → 103 taken 279 times.
|
279 | end subroutine l2_projection_mform_scalar |
| 420 | |||
| 421 | !> @brief Compute the L2 projection of a user-defined function onto a vector-valued m-form space | ||
| 422 | !> | ||
| 423 | !> @param[out] c The resulting m-form containing the L2 projection | ||
| 424 | !> @param[in] space The m-form space onto which the projection is computed | ||
| 425 | !> @param[in] userfun_x The x-component of the user-defined function to project onto the m-form space | ||
| 426 | !> @param[in] userfun_y The y-component of the user-defined function to project onto the m-form space | ||
| 427 | !> @param[in] userfun_z The z-component of the user-defined function to project onto the m-form space | ||
| 428 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 429 | !> @param[in] solver _(optional)_ The solver to use for the L2 projection, if not provided, a solver is created internally | ||
| 430 | !> @param[in] n_quad_extra _(optional)_ The extra number of quadrature points to use for the L2 projection on top of those needed | ||
| 431 | !> to exactly integrate the B-spline spaces | ||
| 432 | !> @param[in] user_fun_is_physical _(optional)_ Whether the user functions are defined in physical coordinates (default: .false.) | ||
| 433 |
1/2✓ Branch 4 → 5 taken 342 times.
✗ Branch 4 → 6 not taken.
|
342 | subroutine l2_projection_mform_vector(c, space, userfun_x, userfun_y, userfun_z, coord_transform, solver, n_quad_extra, & |
| 434 | user_fun_is_physical) | ||
| 435 | use m_tensorprod_linalg, only: l2_projection | ||
| 436 | use m_mform_solver, only: GenericSolver | ||
| 437 | use m_mform_matrix, only: MassMatrix | ||
| 438 | use m_coord_transform_abstract, only: CoordTransformAbstract | ||
| 439 | implicit none | ||
| 440 | |||
| 441 | type(MFormFun), intent(out) :: c | ||
| 442 | type(MFormSpace), intent(in) :: space | ||
| 443 | procedure(user_function_3d_interface) :: userfun_x, userfun_y, userfun_z | ||
| 444 | class(CoordTransformAbstract), intent(in), optional :: coord_transform | ||
| 445 | type(GenericSolver), intent(in), optional :: solver | ||
| 446 | integer, intent(in), optional :: n_quad_extra | ||
| 447 | logical, intent(in), optional :: user_fun_is_physical | ||
| 448 | |||
| 449 | 342 | type(MFormFun) :: f | |
| 450 | 342 | type(MassMatrix) :: mass | |
| 451 |
1/2✓ Branch 2 → 3 taken 342 times.
✗ Branch 2 → 4 not taken.
|
342 | type(GenericSolver) :: solver_ |
| 452 | |||
| 453 |
1/2✗ Branch 6 → 7 not taken.
✓ Branch 6 → 8 taken 342 times.
|
342 | if (space%m == 0 .or. space%m == 3) then |
| 454 | ✗ | error stop "l2_projection for scalar-valued forms requires 1 userfunction" | |
| 455 | end if | ||
| 456 | |||
| 457 | call inner_product(f, space, userfun_x, userfun_y, userfun_z, n_quad_extra=n_quad_extra, & | ||
| 458 | 342 | & user_fun_is_physical=user_fun_is_physical, coord_transform=coord_transform) | |
| 459 | |||
| 460 |
1/2✓ Branch 9 → 10 taken 342 times.
✗ Branch 9 → 12 not taken.
|
342 | if (present(solver)) then |
| 461 | 342 | call solver%solve(c, f) | |
| 462 | else | ||
| 463 | ✗ | call mass%init(space, n_quad_extra=n_quad_extra, coord_transform=coord_transform) | |
| 464 | ✗ | call solver_%init(mass, coord_transform=coord_transform) | |
| 465 | ✗ | call solver_%solve(c, f) | |
| 466 | |||
| 467 | ✗ | call mass%destroy() | |
| 468 | ✗ | call solver_%destroy() | |
| 469 | end if | ||
| 470 | |||
| 471 | 342 | call f%destroy() | |
| 472 | |||
| 473 |
7/80✗ Branch 19 → 20 not taken.
✓ Branch 19 → 21 taken 342 times.
✗ Branch 21 → 22 not taken.
✓ Branch 21 → 23 taken 342 times.
✗ Branch 23 → 24 not taken.
✓ Branch 23 → 95 taken 342 times.
✗ Branch 24 → 25 not taken.
✗ Branch 24 → 26 not taken.
✗ Branch 27 → 28 not taken.
✗ Branch 27 → 71 not taken.
✗ Branch 28 → 29 not taken.
✗ Branch 28 → 31 not taken.
✗ Branch 29 → 30 not taken.
✗ Branch 29 → 31 not taken.
✗ Branch 31 → 32 not taken.
✗ Branch 31 → 70 not taken.
✗ Branch 32 → 33 not taken.
✗ Branch 32 → 34 not taken.
✗ Branch 34 → 35 not taken.
✗ Branch 34 → 47 not taken.
✗ Branch 36 → 37 not taken.
✗ Branch 36 → 46 not taken.
✗ Branch 37 → 38 not taken.
✗ Branch 37 → 39 not taken.
✗ Branch 39 → 40 not taken.
✗ Branch 39 → 41 not taken.
✗ Branch 41 → 42 not taken.
✗ Branch 41 → 43 not taken.
✗ Branch 43 → 44 not taken.
✗ Branch 43 → 45 not taken.
✗ Branch 47 → 48 not taken.
✗ Branch 47 → 58 not taken.
✗ Branch 49 → 50 not taken.
✗ Branch 49 → 57 not taken.
✗ Branch 50 → 51 not taken.
✗ Branch 50 → 52 not taken.
✗ Branch 52 → 53 not taken.
✗ Branch 52 → 54 not taken.
✗ Branch 54 → 55 not taken.
✗ Branch 54 → 56 not taken.
✗ Branch 58 → 59 not taken.
✗ Branch 58 → 69 not taken.
✗ Branch 60 → 61 not taken.
✗ Branch 60 → 68 not taken.
✗ Branch 61 → 62 not taken.
✗ Branch 61 → 63 not taken.
✗ Branch 63 → 64 not taken.
✗ Branch 63 → 65 not taken.
✗ Branch 65 → 66 not taken.
✗ Branch 65 → 67 not taken.
✗ Branch 71 → 72 not taken.
✗ Branch 71 → 73 not taken.
✗ Branch 73 → 74 not taken.
✗ Branch 73 → 75 not taken.
✗ Branch 75 → 76 not taken.
✗ Branch 75 → 77 not taken.
✗ Branch 77 → 78 not taken.
✗ Branch 77 → 79 not taken.
✗ Branch 79 → 80 not taken.
✗ Branch 79 → 81 not taken.
✗ Branch 81 → 82 not taken.
✗ Branch 81 → 94 not taken.
✗ Branch 83 → 84 not taken.
✗ Branch 83 → 93 not taken.
✗ Branch 84 → 85 not taken.
✗ Branch 84 → 86 not taken.
✗ Branch 86 → 87 not taken.
✗ Branch 86 → 88 not taken.
✗ Branch 88 → 89 not taken.
✗ Branch 88 → 90 not taken.
✗ Branch 90 → 91 not taken.
✗ Branch 90 → 92 not taken.
✗ Branch 95 → 96 not taken.
✓ Branch 95 → 97 taken 342 times.
✗ Branch 97 → 98 not taken.
✓ Branch 97 → 99 taken 342 times.
✓ Branch 99 → 100 taken 342 times.
✗ Branch 99 → 101 not taken.
✗ Branch 101 → 102 not taken.
✓ Branch 101 → 103 taken 342 times.
|
342 | end subroutine l2_projection_mform_vector |
| 474 | |||
| 475 | !> @brief Compute the L2 norm of the difference between a scalar-valued m-form and a user-defined function | ||
| 476 | !> | ||
| 477 | !> @param[in] v The m-form to compare with the user-defined function | ||
| 478 | !> @param[in] userfun The user-defined function to compare with the m-form | ||
| 479 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 480 | !> @param[in] n_quad_extra _(optional)_ The extra number of quadrature points to use for the L2 error on top of those needed | ||
| 481 | !> to exactly integrate the B-spline spaces | ||
| 482 | !> @param[in] user_fun_is_physical _(optional)_ Whether the user function is defined in physical coordinates (default: .false.) | ||
| 483 | !> | ||
| 484 | !> @return The L2 norm of the difference between the m-form and the user-defined function | ||
| 485 | 505 | real(wp) function l2_error_mform_scalar(v, userfun, coord_transform, n_quad_extra, user_fun_is_physical) result(ans) | |
| 486 | use m_tensorprod_basis, only: evaluate | ||
| 487 | use m_tensorprod_linalg, only: integrate | ||
| 488 | use m_coord_transform, only: CoordTransformAbstract, jacobian, transform | ||
| 489 | implicit none | ||
| 490 | |||
| 491 | type(MFormFun), intent(in) :: v | ||
| 492 | procedure(user_function_3d_interface) :: userfun | ||
| 493 | class(CoordTransformAbstract), intent(in), optional :: coord_transform | ||
| 494 | integer, intent(in), optional :: n_quad_extra | ||
| 495 | logical, intent(in), optional :: user_fun_is_physical | ||
| 496 | |||
| 497 | logical :: user_fun_is_physical_ | ||
| 498 | |||
| 499 |
1/2✗ Branch 2 → 3 not taken.
✓ Branch 2 → 4 taken 505 times.
|
505 | if (v%space%m == 1 .or. v%space%m == 2) then |
| 500 | ✗ | error stop "l2_error for vector-valued forms requires three userfunctions" | |
| 501 | end if | ||
| 502 | |||
| 503 |
3/4✓ Branch 4 → 5 taken 147 times.
✓ Branch 4 → 7 taken 358 times.
✓ Branch 5 → 6 taken 147 times.
✗ Branch 5 → 7 not taken.
|
505 | if (present(coord_transform)) then |
| 504 | 147 | user_fun_is_physical_ = present(user_fun_is_physical) | |
| 505 | else | ||
| 506 | 358 | user_fun_is_physical_ = .false. | |
| 507 | end if | ||
| 508 | |||
| 509 | 505 | ans = sqrt(integrate(v%space%tp_spaces(1), int_function, n_quad_extra=n_quad_extra)) | |
| 510 | |||
| 511 | contains | ||
| 512 | 51824200 | pure real(wp) function int_function(xp, yp, zp) result(val) | |
| 513 | real(wp), intent(in) :: xp, yp, zp | ||
| 514 | |||
| 515 | real(wp) :: x, y, z, user_val, jac_val | ||
| 516 | |||
| 517 |
3/4✓ Branch 2 → 3 taken 36070744 times.
✓ Branch 2 → 5 taken 15753456 times.
✓ Branch 3 → 4 taken 36070744 times.
✗ Branch 3 → 5 not taken.
|
51824200 | if (present(coord_transform)) then |
| 518 | 36070744 | jac_val = jacobian(coord_transform, xp, yp, zp) | |
| 519 | else | ||
| 520 | jac_val = 1.0_wp | ||
| 521 | end if | ||
| 522 | |||
| 523 |
1/2✓ Branch 5 → 6 taken 51824200 times.
✗ Branch 5 → 7 not taken.
|
51824200 | if (.not. user_fun_is_physical_) then |
| 524 | 51824200 | user_val = userfun(xp, yp, zp) | |
| 525 | else | ||
| 526 | ✗ | if (present(coord_transform)) then | |
| 527 | ✗ | x = transform(coord_transform, 1, xp, yp, zp) | |
| 528 | ✗ | y = transform(coord_transform, 2, xp, yp, zp) | |
| 529 | ✗ | z = transform(coord_transform, 3, xp, yp, zp) | |
| 530 | ✗ | user_val = userfun(x, y, z) | |
| 531 | |||
| 532 | ✗ | if (v%space%m == 3) then | |
| 533 | ! The pullback transformation for 3-forms yields a factor jacobian | ||
| 534 | ✗ | user_val = user_val * jac_val | |
| 535 | end if | ||
| 536 | end if | ||
| 537 | end if | ||
| 538 | 51824200 | val = (user_val - evaluate(v%tp_funs(1), xp, yp, zp))**2 ! value of the logical m-form | |
| 539 | |||
| 540 |
3/4✓ Branch 12 → 13 taken 36070744 times.
✓ Branch 12 → 17 taken 15753456 times.
✓ Branch 13 → 14 taken 36070744 times.
✗ Branch 13 → 17 not taken.
|
51824200 | if (present(coord_transform)) then |
| 541 |
2/2✓ Branch 14 → 15 taken 6758580 times.
✓ Branch 14 → 16 taken 29312164 times.
|
36070744 | if (v%space%m == 3) then |
| 542 | 6758580 | val = val / jac_val | |
| 543 | else | ||
| 544 | 29312164 | val = val * jac_val | |
| 545 | end if | ||
| 546 | end if | ||
| 547 | 51824200 | end function int_function | |
| 548 | end function | ||
| 549 | |||
| 550 | !> @brief Compute the L2 norm of the difference between a vector-valued m-form and a user-defined function | ||
| 551 | !> | ||
| 552 | !> @param[in] v The m-form to compare with the user-defined function | ||
| 553 | !> @param[in] userfun_x The x-component of the user-defined function to compare with the m-form | ||
| 554 | !> @param[in] userfun_y The y-component of the user-defined function to compare with the m-form | ||
| 555 | !> @param[in] userfun_z The z-component of the user-defined function to compare with the m-form | ||
| 556 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 557 | !> @param[in] n_quad_extra _(optional)_ The extra number of quadrature points to use for the L2 error on top of those needed | ||
| 558 | !> to exactly integrate the B-spline spaces | ||
| 559 | !> @param[in] user_fun_is_physical _(optional)_ Whether the user functions are defined in physical coordinates (default: .false.) | ||
| 560 | !> | ||
| 561 | !> @return The L2 norm of the difference between the m-form and the user-defined function | ||
| 562 | 499 | real(wp) function l2_error_mform_vector(v, userfun_x, userfun_y, userfun_z, coord_transform & | |
| 563 | , n_quad_extra, user_fun_is_physical) result(ans) | ||
| 564 | use m_mform_basis, only: MFormSpace, get_zero_form_space | ||
| 565 | use m_tensorprod_basis, only: TensorProdSpace, evaluate | ||
| 566 | use m_tensorprod_linalg, only: integrate | ||
| 567 | use m_coord_transform, only: CoordTransformAbstract, jacobian, transform, jacobian_matrix, jacobian_matrix_inv, G_matrix, & | ||
| 568 | G_matrix_inv | ||
| 569 | implicit none | ||
| 570 | |||
| 571 | type(MFormFun), intent(in) :: v | ||
| 572 | procedure(user_function_3d_interface) :: userfun_x, userfun_y, userfun_z | ||
| 573 | class(CoordTransformAbstract), intent(in), optional :: coord_transform | ||
| 574 | integer, intent(in), optional :: n_quad_extra | ||
| 575 | logical, intent(in), optional :: user_fun_is_physical | ||
| 576 | |||
| 577 | logical :: user_fun_is_physical_ | ||
| 578 | type(MFormSpace) :: space0 | ||
| 579 | |||
| 580 |
1/2✗ Branch 2 → 3 not taken.
✓ Branch 2 → 4 taken 499 times.
|
499 | if (v%space%m == 0 .or. v%space%m == 3) then |
| 581 | ✗ | error stop "l2_error for scalar-valued forms requires one userfunction" | |
| 582 | end if | ||
| 583 | |||
| 584 |
3/4✓ Branch 4 → 5 taken 110 times.
✓ Branch 4 → 7 taken 389 times.
✓ Branch 5 → 6 taken 110 times.
✗ Branch 5 → 7 not taken.
|
499 | if (present(coord_transform)) then |
| 585 | 110 | user_fun_is_physical_ = present(user_fun_is_physical) | |
| 586 | else | ||
| 587 | 389 | user_fun_is_physical_ = .false. | |
| 588 | end if | ||
| 589 | |||
| 590 | ! TODO: user_fun_is_physical (assumed .false. in current implementation) | ||
| 591 | |||
| 592 | 499 | space0 = get_zero_form_space(v%space) | |
| 593 | |||
| 594 |
1/2✓ Branch 10 → 11 taken 499 times.
✗ Branch 10 → 12 not taken.
|
499 | ans = sqrt(integrate(space0%tp_spaces(1), int_function, n_quad_extra=n_quad_extra)) |
| 595 | |||
| 596 | contains | ||
| 597 | 53036460 | pure real(wp) function int_function(xp, yp, zp) result(val) | |
| 598 | real(wp), intent(in) :: xp, yp, zp | ||
| 599 | |||
| 600 | real(wp) :: x, y, z, jac_val, diff_x, diff_y, diff_z | ||
| 601 | real(wp) :: userfun_xp_val, userfun_yp_val, userfun_zp_val, userfun_x_val, userfun_y_val, userfun_z_val | ||
| 602 | real(wp) :: mat(3, 3) | ||
| 603 | |||
| 604 |
1/2✗ Branch 2 → 3 not taken.
✓ Branch 2 → 11 taken 53036460 times.
|
53036460 | if (user_fun_is_physical_) then |
| 605 | ✗ | if (present(coord_transform)) then | |
| 606 | ✗ | x = transform(coord_transform, 1, xp, yp, zp) | |
| 607 | ✗ | y = transform(coord_transform, 2, xp, yp, zp) | |
| 608 | ✗ | z = transform(coord_transform, 3, xp, yp, zp) | |
| 609 | ✗ | userfun_x_val = userfun_x(x, y, z) | |
| 610 | ✗ | userfun_y_val = userfun_y(x, y, z) | |
| 611 | ✗ | userfun_z_val = userfun_z(x, y, z) | |
| 612 | |||
| 613 | ✗ | if (v%space%m == 1) then | |
| 614 | userfun_xp_val = jacobian_matrix(coord_transform, 1, 1, xp, yp, zp) * userfun_x_val + & | ||
| 615 | jacobian_matrix(coord_transform, 2, 1, xp, yp, zp) * userfun_y_val + & | ||
| 616 | ✗ | jacobian_matrix(coord_transform, 3, 1, xp, yp, zp) * userfun_z_val | |
| 617 | userfun_yp_val = jacobian_matrix(coord_transform, 1, 2, xp, yp, zp) * userfun_x_val + & | ||
| 618 | jacobian_matrix(coord_transform, 2, 2, xp, yp, zp) * userfun_y_val + & | ||
| 619 | ✗ | jacobian_matrix(coord_transform, 3, 2, xp, yp, zp) * userfun_z_val | |
| 620 | userfun_zp_val = jacobian_matrix(coord_transform, 1, 3, xp, yp, zp) * userfun_x_val + & | ||
| 621 | jacobian_matrix(coord_transform, 2, 3, xp, yp, zp) * userfun_y_val + & | ||
| 622 | ✗ | jacobian_matrix(coord_transform, 3, 3, xp, yp, zp) * userfun_z_val | |
| 623 | else | ||
| 624 | ✗ | jac_val = jacobian(coord_transform, xp, yp, zp) | |
| 625 | userfun_xp_val = jacobian_matrix_inv(coord_transform, 1, 1, xp, yp, zp) * userfun_x_val + & | ||
| 626 | jacobian_matrix_inv(coord_transform, 1, 2, xp, yp, zp) * userfun_y_val + & | ||
| 627 | ✗ | jacobian_matrix_inv(coord_transform, 1, 3, xp, yp, zp) * userfun_z_val | |
| 628 | userfun_yp_val = jacobian_matrix_inv(coord_transform, 2, 1, xp, yp, zp) * userfun_x_val + & | ||
| 629 | jacobian_matrix_inv(coord_transform, 2, 2, xp, yp, zp) * userfun_y_val + & | ||
| 630 | ✗ | jacobian_matrix_inv(coord_transform, 2, 3, xp, yp, zp) * userfun_z_val | |
| 631 | userfun_zp_val = jacobian_matrix_inv(coord_transform, 3, 1, xp, yp, zp) * userfun_x_val + & | ||
| 632 | jacobian_matrix_inv(coord_transform, 3, 2, xp, yp, zp) * userfun_y_val + & | ||
| 633 | ✗ | jacobian_matrix_inv(coord_transform, 3, 3, xp, yp, zp) * userfun_z_val | |
| 634 | ✗ | userfun_xp_val = userfun_xp_val * jac_val | |
| 635 | ✗ | userfun_yp_val = userfun_yp_val * jac_val | |
| 636 | ✗ | userfun_zp_val = userfun_zp_val * jac_val | |
| 637 | end if | ||
| 638 | end if | ||
| 639 | else | ||
| 640 | 53036460 | userfun_xp_val = userfun_x(xp, yp, zp) | |
| 641 | 53036460 | userfun_yp_val = userfun_y(xp, yp, zp) | |
| 642 | 53036460 | userfun_zp_val = userfun_z(xp, yp, zp) | |
| 643 | end if | ||
| 644 | 53036460 | diff_x = userfun_xp_val - evaluate(v%tp_funs(1), xp, yp, zp) | |
| 645 | 53036460 | diff_y = userfun_yp_val - evaluate(v%tp_funs(2), xp, yp, zp) | |
| 646 | 53036460 | diff_z = userfun_zp_val - evaluate(v%tp_funs(3), xp, yp, zp) | |
| 647 | |||
| 648 |
3/4✓ Branch 14 → 15 taken 26832552 times.
✓ Branch 14 → 16 taken 26203908 times.
✗ Branch 15 → 16 not taken.
✓ Branch 15 → 17 taken 26832552 times.
|
53036460 | if (.not. present(coord_transform)) then |
| 649 | 26203908 | val = diff_x**2 + diff_y**2 + diff_z**2 | |
| 650 | else | ||
| 651 | 26832552 | jac_val = jacobian(coord_transform, xp, yp, zp) | |
| 652 |
2/2✓ Branch 17 → 18 taken 11834580 times.
✓ Branch 17 → 20 taken 14997972 times.
|
26832552 | if (v%space%m == 1) then |
| 653 | 11834580 | mat(1, 1) = G_matrix_inv(coord_transform, 1, 1, xp, yp, zp) * jac_val | |
| 654 | 11834580 | mat(2, 2) = G_matrix_inv(coord_transform, 2, 2, xp, yp, zp) * jac_val | |
| 655 | 11834580 | mat(3, 3) = G_matrix_inv(coord_transform, 3, 3, xp, yp, zp) * jac_val | |
| 656 |
2/2✓ Branch 18 → 19 taken 2644416 times.
✓ Branch 18 → 22 taken 9190164 times.
|
11834580 | if (.not. coord_transform%is_orthogonal) then |
| 657 | 2644416 | mat(1, 2) = G_matrix_inv(coord_transform, 1, 2, xp, yp, zp) * jac_val | |
| 658 | 2644416 | mat(1, 3) = G_matrix_inv(coord_transform, 1, 3, xp, yp, zp) * jac_val | |
| 659 | 2644416 | mat(2, 3) = G_matrix_inv(coord_transform, 2, 3, xp, yp, zp) * jac_val | |
| 660 | end if | ||
| 661 | else ! if (v%space%m == 2) then | ||
| 662 | 14997972 | mat(1, 1) = G_matrix(coord_transform, 1, 1, xp, yp, zp) / jac_val | |
| 663 | 14997972 | mat(2, 2) = G_matrix(coord_transform, 2, 2, xp, yp, zp) / jac_val | |
| 664 | 14997972 | mat(3, 3) = G_matrix(coord_transform, 3, 3, xp, yp, zp) / jac_val | |
| 665 |
2/2✓ Branch 20 → 21 taken 2644416 times.
✓ Branch 20 → 22 taken 12353556 times.
|
14997972 | if (.not. coord_transform%is_orthogonal) then |
| 666 | 2644416 | mat(1, 2) = G_matrix(coord_transform, 1, 2, xp, yp, zp) / jac_val | |
| 667 | 2644416 | mat(1, 3) = G_matrix(coord_transform, 1, 3, xp, yp, zp) / jac_val | |
| 668 | 2644416 | mat(2, 3) = G_matrix(coord_transform, 2, 3, xp, yp, zp) / jac_val | |
| 669 | end if | ||
| 670 | end if | ||
| 671 | 26832552 | val = diff_x**2 * mat(1, 1) + diff_y**2 * mat(2, 2) + diff_z**2 * mat(3, 3) | |
| 672 |
2/2✓ Branch 22 → 23 taken 5288832 times.
✓ Branch 22 → 24 taken 21543720 times.
|
26832552 | if (.not. coord_transform%is_orthogonal) then |
| 673 | 5288832 | val = val + 2.0_wp * (diff_x * diff_y * mat(1, 2) + diff_x * diff_z * mat(1, 3) + diff_y * diff_z * mat(2, 3)) | |
| 674 | end if | ||
| 675 | end if | ||
| 676 | |||
| 677 | 53036460 | end function int_function | |
| 678 | |||
| 679 | end function | ||
| 680 | |||
| 681 | !> @brief Compute the L2 norm of an m-form | ||
| 682 | !> | ||
| 683 | !> @param[in] v The m-form | ||
| 684 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 685 | !> | ||
| 686 | !> @return The L2 norm of the m-form | ||
| 687 | 39 | real(wp) function l2_norm_mform(v, coord_transform) result(ans) | |
| 688 | use m_common, only: zero_function | ||
| 689 | use m_coord_transform, only: CoordTransformAbstract | ||
| 690 | implicit none | ||
| 691 | |||
| 692 | type(MFormFun), intent(in) :: v | ||
| 693 | class(CoordTransformAbstract), intent(in), optional :: coord_transform | ||
| 694 | |||
| 695 |
2/2✓ Branch 2 → 3 taken 37 times.
✓ Branch 2 → 4 taken 2 times.
|
39 | if (v%space%m == 1 .or. v%space%m == 2) then |
| 696 | 37 | ans = l2_error(v, zero_function, zero_function, zero_function, coord_transform=coord_transform) | |
| 697 | else | ||
| 698 | 2 | ans = l2_error(v, zero_function, coord_transform=coord_transform) | |
| 699 | end if | ||
| 700 | 39 | end function | |
| 701 | |||
| 702 | !> @brief Find the weakly divergence-free one-form whose curl is the given two-form (in 2D, assuming homogeneous BC) | ||
| 703 | !> | ||
| 704 | !> @param[out] oneform The resulting one-form whose curl is the given two-form | ||
| 705 | !> @param[in] twoform_x The x-component of the two-form whose curl 'inverse' is to be computed | ||
| 706 | !> @param[in] twoform_y The y-component of the two-form whose curl 'inverse' is to be computed | ||
| 707 | !> @param[in] twoform_space The space of the provided two-form components | ||
| 708 | !> @param[in] cylinder_transform The coordinate transformation associated with the m-form space | ||
| 709 | !> @param[inout] zeroform_solver _(optional)_ The divgrad solver, if not provided, a solver is created internally | ||
| 710 | !> @param[in] constraint _(optional)_ The constraint to be applied to both solvers (besides the boundary conditions) | ||
| 711 |
2/2✓ Branch 4 → 5 taken 235 times.
✓ Branch 4 → 6 taken 84 times.
|
319 | subroutine solve_curlA_equals_B_tensorprod(oneform, twoform_x, twoform_y, twoform_space, cylinder_transform, zeroform_solver, & |
| 712 | constraint) | ||
| 713 | use m_common, only: zero_function, user_function_3d_interface, cumsum | ||
| 714 | use m_bspline_basis, only: BSplineFun | ||
| 715 | use m_mform_solver, only: GenericSolver | ||
| 716 | use m_mform_matrix, only: MassMatrix, DiffDiffMatrix | ||
| 717 | use m_mform_basis, only: MFormFun, gradient, gradient_adjoint, curl_adjoint, previous | ||
| 718 | use m_mform_constraint_abstract, only: MFormConstraintLocal | ||
| 719 | use m_mform_constraint_boundary, only: MFormDirichlet, MFormNeumann | ||
| 720 | use m_coord_transform, only: CylinderTransform | ||
| 721 | implicit none | ||
| 722 | |||
| 723 | type(MFormFun), intent(out) :: oneform | ||
| 724 | type(BSplineFun), intent(in) :: twoform_x, twoform_y | ||
| 725 | type(MFormSpace), intent(in) :: twoform_space | ||
| 726 | type(GenericSolver), intent(inout), target, optional :: zeroform_solver | ||
| 727 | class(MFormConstraintLocal), intent(in), optional :: constraint | ||
| 728 | type(CylinderTransform), intent(in), optional :: cylinder_transform | ||
| 729 | 319 | type(GenericSolver), target :: zeroform_solver_ | |
| 730 | type(GenericSolver), pointer :: zeroform_solver_p | ||
| 731 | 319 | type(DiffDiffMatrix) :: divgrad | |
| 732 |
2/2✓ Branch 2 → 3 taken 238 times.
✓ Branch 2 → 4 taken 81 times.
|
319 | type(MFormFun) :: f, df, zeroform, grad, oneform_tmp |
| 733 | type(MFormSpace) :: oneform_space | ||
| 734 | class(MFormConstraintLocal), allocatable :: constraint_local | ||
| 735 | integer :: twoform_m, i, j, k | ||
| 736 | character(len=100) :: err_msg | ||
| 737 | |||
| 738 | ! TODO calling previous() multiple times (for each twoform) is not efficient | ||
| 739 | 319 | oneform_space = previous(twoform_space) | |
| 740 | 319 | twoform_m = twoform_space%m | |
| 741 | |||
| 742 |
1/2✗ Branch 7 → 8 not taken.
✓ Branch 7 → 14 taken 319 times.
|
319 | if (twoform_m /= 2) then |
| 743 | ✗ | write (err_msg, '(A,I0,A)') "ERROR in solve_curlA_equals_B_ori: the provided space is a ", twoform_m, & | |
| 744 | ✗ | "-form, but a 2-form is required." | |
| 745 | ✗ | error stop err_msg | |
| 746 | end if | ||
| 747 |
1/2✗ Branch 14 → 15 not taken.
✓ Branch 14 → 16 taken 319 times.
|
319 | if (twoform_space%dimensionality == 3) then |
| 748 | ✗ | error stop "ERROR in solve_curlA_equals_B_ori: the provided space is not 2D, but a 2D space is required." | |
| 749 | end if | ||
| 750 | |||
| 751 |
2/2✓ Branch 16 → 17 taken 35 times.
✓ Branch 16 → 60 taken 284 times.
|
319 | if (present(zeroform_solver)) then |
| 752 | zeroform_solver_p => zeroform_solver | ||
| 753 | else | ||
| 754 |
1/2✓ Branch 19 → 20 taken 35 times.
✗ Branch 19 → 21 not taken.
|
35 | call divgrad%init(previous(oneform_space), coord_transform=cylinder_transform) |
| 755 |
4/34✗ Branch 22 → 23 not taken.
✓ Branch 22 → 24 taken 35 times.
✗ Branch 24 → 25 not taken.
✓ Branch 24 → 37 taken 35 times.
✗ Branch 26 → 27 not taken.
✗ Branch 26 → 36 not taken.
✗ Branch 27 → 28 not taken.
✗ Branch 27 → 29 not taken.
✗ Branch 29 → 30 not taken.
✗ Branch 29 → 31 not taken.
✗ Branch 31 → 32 not taken.
✗ Branch 31 → 33 not taken.
✗ Branch 33 → 34 not taken.
✗ Branch 33 → 35 not taken.
✗ Branch 37 → 38 not taken.
✓ Branch 37 → 48 taken 35 times.
✗ Branch 39 → 40 not taken.
✗ Branch 39 → 47 not taken.
✗ Branch 40 → 41 not taken.
✗ Branch 40 → 42 not taken.
✗ Branch 42 → 43 not taken.
✗ Branch 42 → 44 not taken.
✗ Branch 44 → 45 not taken.
✗ Branch 44 → 46 not taken.
✗ Branch 48 → 49 not taken.
✓ Branch 48 → 59 taken 35 times.
✗ Branch 50 → 51 not taken.
✗ Branch 50 → 58 not taken.
✗ Branch 51 → 52 not taken.
✗ Branch 51 → 53 not taken.
✗ Branch 53 → 54 not taken.
✗ Branch 53 → 55 not taken.
✗ Branch 55 → 56 not taken.
✗ Branch 55 → 57 not taken.
|
35 | call zeroform_solver_%init(divgrad, constraint1=MFormNeumann(), constraint2=constraint, coord_transform=cylinder_transform) |
| 756 | zeroform_solver_p => zeroform_solver_ | ||
| 757 | end if | ||
| 758 | |||
| 759 | ! Step one: find a one-form such that curl oneform = twoform | ||
| 760 | ! This can be done analytically since the two-form is given as a product of 1D B-splines | ||
| 761 | 319 | call twoform_x_integrated%init(oneform_space%tp_spaces(2)%spaces(1)) | |
| 762 | 319 | call cumsum(twoform_x_integrated%data, twoform_x%data) | |
| 763 | 319 | call twoform_y_copy%init(twoform_y%bspline) | |
| 764 |
4/10✓ Branch 63 → 64 taken 319 times.
✗ Branch 63 → 65 not taken.
✗ Branch 64 → 65 not taken.
✓ Branch 64 → 71 taken 319 times.
✗ Branch 65 → 66 not taken.
✗ Branch 65 → 67 not taken.
✗ Branch 67 → 68 not taken.
✗ Branch 67 → 69 not taken.
✓ Branch 72 → 73 taken 10057 times.
✓ Branch 72 → 74 taken 319 times.
|
10695 | twoform_y_copy%data = twoform_y%data |
| 765 | 319 | call oneform_tmp%init(oneform_space) | |
| 766 |
2/2✓ Branch 76 → 77 taken 319 times.
✓ Branch 76 → 86 taken 319 times.
|
638 | do k = oneform_space%tp_spaces(2)%rank_resp_bspline%k0, oneform_space%tp_spaces(2)%rank_resp_bspline%k1 |
| 767 |
2/2✓ Branch 78 → 79 taken 9116 times.
✓ Branch 78 → 84 taken 319 times.
|
9754 | do j = oneform_space%tp_spaces(2)%rank_resp_bspline%j0, oneform_space%tp_spaces(2)%rank_resp_bspline%j1 |
| 768 |
2/2✓ Branch 80 → 81 taken 212668 times.
✓ Branch 80 → 82 taken 9116 times.
|
222103 | do i = oneform_space%tp_spaces(2)%rank_resp_bspline%i0, oneform_space%tp_spaces(2)%rank_resp_bspline%i1 |
| 769 | 221784 | oneform_tmp%tp_funs(2)%data(i, j, k) = twoform_x_integrated%data(i) * twoform_y%data(j) | |
| 770 | end do | ||
| 771 | end do | ||
| 772 | end do | ||
| 773 | |||
| 774 | 319 | call oneform_tmp%distribute() | |
| 775 | |||
| 776 | ! Step two: project the one-form onto the weakly divergence-free space | ||
| 777 | 319 | call inner_product(f, oneform_space, zero_function, wrapped_oneform_y, zero_function, cylinder_transform) | |
| 778 | |||
| 779 | 319 | call gradient_adjoint(df, f) | |
| 780 | 319 | call zeroform_solver_p%solve(zeroform, df) | |
| 781 | 319 | call gradient(grad, zeroform) | |
| 782 | 319 | call oneform_tmp%axpy(-1.0_wp, grad) | |
| 783 | |||
| 784 |
3/4✓ Branch 93 → 94 taken 35 times.
✓ Branch 93 → 107 taken 284 times.
✓ Branch 94 → 95 taken 35 times.
✗ Branch 94 → 107 not taken.
|
319 | if (present(constraint)) then |
| 785 |
1/2✗ Branch 95 → 96 not taken.
✓ Branch 95 → 97 taken 35 times.
|
35 | allocate (constraint_local, source=constraint) |
| 786 | 35 | call constraint_local%init(oneform_space, cylinder_transform) | |
| 787 | 35 | call constraint_local%apply(oneform, oneform_tmp) | |
| 788 | 35 | call constraint_local%destroy() | |
| 789 |
2/4✗ Branch 101 → 102 not taken.
✓ Branch 101 → 103 taken 35 times.
✓ Branch 103 → 104 taken 35 times.
✗ Branch 103 → 106 not taken.
|
70 | deallocate (constraint_local) |
| 790 | else | ||
| 791 | ! NOTE: this is a deep copy using the overloaded assignment(=) operator | ||
| 792 | 284 | oneform = oneform_tmp | |
| 793 | end if | ||
| 794 | |||
| 795 | 319 | call f%destroy() | |
| 796 | 319 | call df%destroy() | |
| 797 | 319 | call grad%destroy() | |
| 798 | 319 | call zeroform%destroy() | |
| 799 | 319 | call oneform_tmp%destroy() | |
| 800 | 319 | call twoform_x_integrated%destroy() | |
| 801 | |||
| 802 |
21/98✓ Branch 115 → 116 taken 35 times.
✓ Branch 115 → 119 taken 284 times.
✓ Branch 119 → 120 taken 35 times.
✓ Branch 119 → 121 taken 284 times.
✓ Branch 121 → 122 taken 35 times.
✓ Branch 121 → 123 taken 284 times.
✗ Branch 123 → 124 not taken.
✓ Branch 123 → 195 taken 319 times.
✗ Branch 124 → 125 not taken.
✗ Branch 124 → 126 not taken.
✗ Branch 127 → 128 not taken.
✗ Branch 127 → 171 not taken.
✗ Branch 128 → 129 not taken.
✗ Branch 128 → 131 not taken.
✗ Branch 129 → 130 not taken.
✗ Branch 129 → 131 not taken.
✗ Branch 131 → 132 not taken.
✗ Branch 131 → 170 not taken.
✗ Branch 132 → 133 not taken.
✗ Branch 132 → 134 not taken.
✗ Branch 134 → 135 not taken.
✗ Branch 134 → 147 not taken.
✗ Branch 136 → 137 not taken.
✗ Branch 136 → 146 not taken.
✗ Branch 137 → 138 not taken.
✗ Branch 137 → 139 not taken.
✗ Branch 139 → 140 not taken.
✗ Branch 139 → 141 not taken.
✗ Branch 141 → 142 not taken.
✗ Branch 141 → 143 not taken.
✗ Branch 143 → 144 not taken.
✗ Branch 143 → 145 not taken.
✗ Branch 147 → 148 not taken.
✗ Branch 147 → 158 not taken.
✗ Branch 149 → 150 not taken.
✗ Branch 149 → 157 not taken.
✗ Branch 150 → 151 not taken.
✗ Branch 150 → 152 not taken.
✗ Branch 152 → 153 not taken.
✗ Branch 152 → 154 not taken.
✗ Branch 154 → 155 not taken.
✗ Branch 154 → 156 not taken.
✗ Branch 158 → 159 not taken.
✗ Branch 158 → 169 not taken.
✗ Branch 160 → 161 not taken.
✗ Branch 160 → 168 not taken.
✗ Branch 161 → 162 not taken.
✗ Branch 161 → 163 not taken.
✗ Branch 163 → 164 not taken.
✗ Branch 163 → 165 not taken.
✗ Branch 165 → 166 not taken.
✗ Branch 165 → 167 not taken.
✗ Branch 171 → 172 not taken.
✗ Branch 171 → 173 not taken.
✗ Branch 173 → 174 not taken.
✗ Branch 173 → 175 not taken.
✗ Branch 175 → 176 not taken.
✗ Branch 175 → 177 not taken.
✗ Branch 177 → 178 not taken.
✗ Branch 177 → 179 not taken.
✗ Branch 179 → 180 not taken.
✗ Branch 179 → 181 not taken.
✗ Branch 181 → 182 not taken.
✗ Branch 181 → 194 not taken.
✗ Branch 183 → 184 not taken.
✗ Branch 183 → 193 not taken.
✗ Branch 184 → 185 not taken.
✗ Branch 184 → 186 not taken.
✗ Branch 186 → 187 not taken.
✗ Branch 186 → 188 not taken.
✗ Branch 188 → 189 not taken.
✗ Branch 188 → 190 not taken.
✗ Branch 190 → 191 not taken.
✗ Branch 190 → 192 not taken.
✓ Branch 195 → 196 taken 319 times.
✗ Branch 195 → 197 not taken.
✗ Branch 197 → 198 not taken.
✓ Branch 197 → 199 taken 319 times.
✓ Branch 199 → 200 taken 319 times.
✗ Branch 199 → 201 not taken.
✗ Branch 201 → 202 not taken.
✓ Branch 201 → 203 taken 319 times.
✓ Branch 203 → 204 taken 319 times.
✗ Branch 203 → 205 not taken.
✗ Branch 205 → 206 not taken.
✓ Branch 205 → 207 taken 319 times.
✓ Branch 207 → 208 taken 319 times.
✗ Branch 207 → 209 not taken.
✗ Branch 209 → 210 not taken.
✓ Branch 209 → 211 taken 319 times.
✓ Branch 211 → 212 taken 35 times.
✓ Branch 211 → 213 taken 284 times.
✓ Branch 213 → 214 taken 35 times.
✓ Branch 213 → 215 taken 284 times.
✓ Branch 215 → 216 taken 319 times.
✗ Branch 215 → 217 not taken.
✗ Branch 217 → 218 not taken.
✓ Branch 217 → 219 taken 319 times.
|
354 | if (.not. present(zeroform_solver)) then |
| 803 | 35 | call divgrad%destroy() | |
| 804 | 35 | call zeroform_solver_%destroy() | |
| 805 | end if | ||
| 806 | contains | ||
| 807 | ! pure real(wp) function wrapped_oneform_x(xp, yp, zp) result(ans) | ||
| 808 | ! implicit none | ||
| 809 | |||
| 810 | ! real(wp), intent(in) :: xp, yp, zp | ||
| 811 | |||
| 812 | ! ans = evaluate(oneform_tmp, 1, xp, yp, zp, eval_logical_form=.true.) | ||
| 813 | ! end function wrapped_oneform_x | ||
| 814 | |||
| 815 | 23375904 | pure real(wp) function wrapped_oneform_y(xp, yp, zp) result(ans) | |
| 816 | use m_bspline_basis, only: evaluate | ||
| 817 | implicit none | ||
| 818 | |||
| 819 | real(wp), intent(in) :: xp, yp, zp | ||
| 820 | |||
| 821 | ! ans = evaluate(oneform_tmp, 2, xp, yp, zp, eval_logical_form=.true.) | ||
| 822 | 23375904 | ans = evaluate(twoform_x_integrated, xp) * evaluate(twoform_y_copy, yp) | |
| 823 | 23375904 | end function wrapped_oneform_y | |
| 824 | end subroutine solve_curlA_equals_B_tensorprod | ||
| 825 | |||
| 826 | end module m_mform_linalg | ||
| 827 |