mform/m_mform_matrix.F90
| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | !> @brief Module for m-form matrices | ||
| 2 | !> | ||
| 3 | !> This module provides mass matrices for each of the m-form spaces, as well as matrix representations of the diffrential operators | ||
| 4 | module m_mform_matrix | ||
| 5 | use m_mform_basis, only: MFormSpace | ||
| 6 | #include "petsc.fi" | ||
| 7 | #include "timer.fi" | ||
| 8 | |||
| 9 | implicit none | ||
| 10 | |||
| 11 | private | ||
| 12 | public :: AbstractMatrix, MassMatrix, DiffMatrix, DiffDiffMatrix, ZeroMatrix | ||
| 13 | public :: gradient_matrix, curl_matrix, divergence_matrix, transpose | ||
| 14 | |||
| 15 | !> @brief Abstract base type for m-form matrices | ||
| 16 | type, abstract :: AbstractMatrix | ||
| 17 | !> The m-form space of the solution vector | ||
| 18 | type(MFormSpace) :: space_col | ||
| 19 | |||
| 20 | !> The m-form space of the residual | ||
| 21 | type(MFormSpace) :: space_row | ||
| 22 | |||
| 23 | !> PETSc matrix | ||
| 24 | Mat :: mat | ||
| 25 | |||
| 26 | !> Whether the matrix is symmetric | ||
| 27 | logical :: is_symmetric | ||
| 28 | !> Whether the matrix is symmetric positive definite | ||
| 29 | logical :: is_spd | ||
| 30 | contains | ||
| 31 | procedure(constraint_coeff_interface), deferred :: constraint_coeff | ||
| 32 | procedure :: destroy => destroy_abstract_matrix | ||
| 33 | end type | ||
| 34 | |||
| 35 | interface | ||
| 36 | pure function constraint_coeff_interface(this) result(ans) | ||
| 37 | use m_common, only: wp | ||
| 38 | import :: AbstractMatrix | ||
| 39 | class(AbstractMatrix), intent(in) :: this | ||
| 40 | |||
| 41 | real(wp) :: ans | ||
| 42 | end function constraint_coeff_interface | ||
| 43 | end interface | ||
| 44 | |||
| 45 | !> @brief Mass matrix | ||
| 46 | type, extends(AbstractMatrix) :: MassMatrix | ||
| 47 | |||
| 48 | contains | ||
| 49 | procedure :: init => init_mass_matrix | ||
| 50 | procedure :: constraint_coeff => constraint_coeff_mass | ||
| 51 | procedure :: matvec => matvec_mass | ||
| 52 | end type | ||
| 53 | |||
| 54 | !> @brief Matrix representation of the weak form of a differential operator (gradient, curl, divergence) | ||
| 55 | !> | ||
| 56 | !> The matrix is given by M * D, where M is the mass matrix of the next m-form space and D is the 'finite difference' matrix | ||
| 57 | type, extends(AbstractMatrix) :: DiffMatrix | ||
| 58 | contains | ||
| 59 | procedure :: init => init_diff_matrix | ||
| 60 | procedure :: constraint_coeff => constraint_coeff_diff | ||
| 61 | end type | ||
| 62 | |||
| 63 | !> @brief Matrix representation of the weak form of a differential operator (divgrad, curlcurl, graddiv) | ||
| 64 | !> | ||
| 65 | !> The matrix is given by D^T * M * D, where M is the mass matrix of the next m-form space and D is the 'finite difference' matrix | ||
| 66 | type, extends(AbstractMatrix) :: DiffDiffMatrix | ||
| 67 | |||
| 68 | contains | ||
| 69 | procedure :: init => init_diffdiff_matrix | ||
| 70 | procedure :: constraint_coeff => constraint_coeff_diffdiff | ||
| 71 | end type | ||
| 72 | |||
| 73 | type, extends(AbstractMatrix) :: ZeroMatrix | ||
| 74 | contains | ||
| 75 | procedure :: init => init_zero_matrix | ||
| 76 | procedure :: constraint_coeff => constraint_coeff_zero | ||
| 77 | end type | ||
| 78 | |||
| 79 | interface transpose | ||
| 80 | module procedure transpose_abstract_matrix | ||
| 81 | end interface transpose | ||
| 82 | |||
| 83 | contains | ||
| 84 | |||
| 85 | !> @brief Initialize the mass matrix for the given m-form space | ||
| 86 | !> | ||
| 87 | !> @param[inout] this The MassMatrix object to initialize | ||
| 88 | !> @param[in] space The MFormSpace for which to initialize the mass matrix | ||
| 89 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 90 | !> @param[in] userfun _(optional)_ A user-defined function to used to scale the inner product that defines the mass matrix | ||
| 91 | !> @param[in] userfun_mat _(optional)_ A user-defined matrix-valued function to used to scale the inner product that defines the | ||
| 92 | !> mass matrix | ||
| 93 | !> @param[in] n_quad_extra _(optional)_ The number of extra quadrature points to use for the mass matrix on top of those needed | ||
| 94 | !> for the exact integration of the product of the B-spline basis functions | ||
| 95 | 613 | subroutine init_mass_matrix(this, space, coord_transform, userfun, userfun_mat, n_quad_extra) | |
| 96 | use m_common, only: user_function_3d_interface, user_function_3d_mat_interface | ||
| 97 | use m_coord_transform, only: CoordTransformAbstract | ||
| 98 | |||
| 99 | implicit none | ||
| 100 | class(MassMatrix), intent(inout) :: this | ||
| 101 | type(MFormSpace), intent(in) :: space | ||
| 102 | class(CoordTransformAbstract), optional, intent(in) :: coord_transform | ||
| 103 | procedure(user_function_3d_interface), optional :: userfun | ||
| 104 | procedure(user_function_3d_mat_interface), optional :: userfun_mat | ||
| 105 | integer, intent(in), optional :: n_quad_extra | ||
| 106 | |||
| 107 | 613 | call this%destroy() | |
| 108 | |||
| 109 |
5/8✓ Branch 3 → 4 taken 613 times.
✗ Branch 3 → 7 not taken.
✓ Branch 4 → 5 taken 613 times.
✗ Branch 4 → 6 not taken.
✓ Branch 7 → 8 taken 613 times.
✗ Branch 7 → 10 not taken.
✓ Branch 8 → 9 taken 309 times.
✓ Branch 8 → 10 taken 304 times.
|
613 | this%space_col = space |
| 110 |
5/8✓ Branch 10 → 11 taken 613 times.
✗ Branch 10 → 14 not taken.
✓ Branch 11 → 12 taken 613 times.
✗ Branch 11 → 13 not taken.
✓ Branch 14 → 15 taken 613 times.
✗ Branch 14 → 17 not taken.
✓ Branch 15 → 16 taken 309 times.
✓ Branch 15 → 17 taken 304 times.
|
613 | this%space_row = space |
| 111 | |||
| 112 |
2/2✓ Branch 17 → 18 taken 155 times.
✓ Branch 17 → 19 taken 458 times.
|
613 | if (space%m == 0 .or. space%m == 3) then |
| 113 | 155 | call init_mass_matrix_scalar(this, space, coord_transform=coord_transform, userfun=userfun, n_quad_extra=n_quad_extra) | |
| 114 | else | ||
| 115 | call init_mass_matrix_vector(this, space, coord_transform=coord_transform, userfun=userfun, userfun_mat=userfun_mat, & | ||
| 116 | 458 | n_quad_extra=n_quad_extra) | |
| 117 | end if | ||
| 118 | |||
| 119 | 613 | this%is_symmetric = .true. | |
| 120 | 613 | this%is_spd = .true. | |
| 121 | 613 | end subroutine | |
| 122 | |||
| 123 | !> @brief Initialize the mass matrix for a scalar m-form space (i.e., m=0,3) | ||
| 124 | !> | ||
| 125 | !> @param[inout] this The MassMatrix object to initialize | ||
| 126 | !> @param[in] space The MFormSpace for which to initialize the mass matrix | ||
| 127 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 128 | !> @param[in] userfun _(optional)_ A user-defined function to used to scale the inner product that defines the mass matrix | ||
| 129 | !> @param[in] n_quad_extra _(optional)_ The number of extra quadrature points to use for the mass matrix on top of those needed | ||
| 130 | !> for the exact integration of the product of the B-spline basis functions | ||
| 131 | 155 | subroutine init_mass_matrix_scalar(this, space, coord_transform, userfun, n_quad_extra) | |
| 132 | use m_tensorprod_matrix, only: TensorProdMat, mass_matrix, get_petsc | ||
| 133 | use m_common, only: user_function_3d_interface | ||
| 134 | use m_coord_transform, only: CoordTransformAbstract, jacobian | ||
| 135 | |||
| 136 | implicit none | ||
| 137 | class(MassMatrix), intent(inout) :: this | ||
| 138 | type(MFormSpace), intent(in) :: space | ||
| 139 | class(CoordTransformAbstract), optional, intent(in) :: coord_transform | ||
| 140 | procedure(user_function_3d_interface), optional :: userfun | ||
| 141 | integer, intent(in), optional :: n_quad_extra | ||
| 142 | |||
| 143 | 155 | type(TensorProdMat) :: M_TP | |
| 144 | |||
| 145 | TIMER_DECL(t_total) | ||
| 146 | TIMER_DECL(t_get_petsc) | ||
| 147 | TIMER_INIT_START(t_total, "MassMatrix::init") | ||
| 148 | |||
| 149 |
4/6✓ Branch 2 → 3 taken 67 times.
✓ Branch 2 → 4 taken 88 times.
✗ Branch 3 → 4 not taken.
✓ Branch 3 → 5 taken 67 times.
✗ Branch 4 → 5 not taken.
✓ Branch 4 → 8 taken 88 times.
|
155 | if (present(coord_transform) .or. present(userfun)) then |
| 150 |
2/2✓ Branch 5 → 6 taken 49 times.
✓ Branch 5 → 7 taken 18 times.
|
67 | if (space%m == 0) then |
| 151 | 49 | call mass_matrix(M_TP, space%tp_spaces(1), space%tp_spaces(1), coordfun=jacobian_wrap, n_quad_extra=n_quad_extra) | |
| 152 | else ! if (space%m == 3) then | ||
| 153 | 18 | call mass_matrix(M_TP, space%tp_spaces(1), space%tp_spaces(1), coordfun=jacobian_inverse_wrap, n_quad_extra=n_quad_extra) | |
| 154 | end if | ||
| 155 | else | ||
| 156 | 88 | call mass_matrix(M_TP, space%tp_spaces(1), space%tp_spaces(1)) | |
| 157 | end if | ||
| 158 | |||
| 159 | TIMER_INIT_START(t_get_petsc, "MassMatrix::get_petsc") | ||
| 160 | 155 | call get_petsc(this%mat, M_TP) | |
| 161 | TIMER_STOP(t_get_petsc) | ||
| 162 | |||
| 163 |
3/6✗ Branch 11 → 12 not taken.
✓ Branch 11 → 13 taken 155 times.
✗ Branch 13 → 14 not taken.
✓ Branch 13 → 15 taken 155 times.
✓ Branch 15 → 16 taken 155 times.
✗ Branch 15 → 17 not taken.
|
155 | call M_TP%destroy() |
| 164 | |||
| 165 | TIMER_STOP(t_total) | ||
| 166 | contains | ||
| 167 | 11965784 | pure real(wp) function jacobian_wrap(xp, yp, zp) | |
| 168 | use m_common, only: wp | ||
| 169 | real(wp), intent(in) :: xp, yp, zp | ||
| 170 | |||
| 171 |
2/4✓ Branch 2 → 3 taken 11965784 times.
✗ Branch 2 → 5 not taken.
✓ Branch 3 → 4 taken 11965784 times.
✗ Branch 3 → 5 not taken.
|
11965784 | if (.not. present(coord_transform)) then |
| 172 | jacobian_wrap = 1._wp | ||
| 173 | else | ||
| 174 | 11965784 | jacobian_wrap = jacobian(coord_transform, xp, yp, zp) | |
| 175 | end if | ||
| 176 | |||
| 177 |
1/2✗ Branch 5 → 6 not taken.
✓ Branch 5 → 8 taken 11965784 times.
|
11965784 | if (present(userfun)) then |
| 178 | ✗ | jacobian_wrap = jacobian_wrap * userfun(xp, yp, zp) | |
| 179 | end if | ||
| 180 | |||
| 181 | 11965784 | end function jacobian_wrap | |
| 182 | 6431964 | pure real(wp) function jacobian_inverse_wrap(xp, yp, zp) | |
| 183 | use m_common, only: wp | ||
| 184 | real(wp), intent(in) :: xp, yp, zp | ||
| 185 | |||
| 186 |
2/4✓ Branch 2 → 3 taken 6431964 times.
✗ Branch 2 → 5 not taken.
✓ Branch 3 → 4 taken 6431964 times.
✗ Branch 3 → 5 not taken.
|
6431964 | if (.not. present(coord_transform)) then |
| 187 | jacobian_inverse_wrap = 1._wp | ||
| 188 | else | ||
| 189 | 6431964 | jacobian_inverse_wrap = 1._wp / jacobian(coord_transform, xp, yp, zp) | |
| 190 | end if | ||
| 191 | |||
| 192 |
1/2✗ Branch 5 → 6 not taken.
✓ Branch 5 → 8 taken 6431964 times.
|
6431964 | if (present(userfun)) then |
| 193 | ✗ | jacobian_inverse_wrap = jacobian_inverse_wrap * userfun(xp, yp, zp) | |
| 194 | end if | ||
| 195 | |||
| 196 | 6431964 | end function jacobian_inverse_wrap | |
| 197 | |||
| 198 | end subroutine | ||
| 199 | |||
| 200 | !> @brief Initialize the mass matrix for a vector-valued m-form space (i.e., m=1,2) | ||
| 201 | !> | ||
| 202 | !> @param[inout] this The MassMatrix object to initialize | ||
| 203 | !> @param[in] space The MFormSpace for which to initialize the mass matrix | ||
| 204 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 205 | !> @param[in] userfun _(optional)_ A user-defined function to used to scale the inner product that defines the mass matrix | ||
| 206 | !> @param[in] userfun_mat _(optional)_ A user-defined matrix-valued function to used to scale the inner product that defines the | ||
| 207 | !> mass matrix | ||
| 208 | !> @param[in] n_quad_extra _(optional)_ The number of extra quadrature points to use for the mass matrix on top of those needed | ||
| 209 | !> for the exact integration of the product of the B-spline basis functions | ||
| 210 | !> | ||
| 211 | !> @note If both userfun and userfun_mat are present, the product of the two is used to scale the inner product defining the mass | ||
| 212 | !> matrix | ||
| 213 | 458 | subroutine init_mass_matrix_vector(this, space, coord_transform, userfun, userfun_mat, n_quad_extra) | |
| 214 | use m_tensorprod_basis, only: size | ||
| 215 | use m_tensorprod_matrix, only: TensorProdMat, mass_matrix, get_petsc | ||
| 216 | use m_common, only: user_function_3d_interface, user_function_3d_mat_interface | ||
| 217 | use m_coord_transform, only: CoordTransformAbstract, jacobian, jacobian_matrix, jacobian_matrix_inv, G_matrix, G_matrix_inv | ||
| 218 | |||
| 219 | implicit none | ||
| 220 | class(MassMatrix), intent(inout) :: this | ||
| 221 | type(MFormSpace), intent(in) :: space | ||
| 222 | class(CoordTransformAbstract), optional, intent(in) :: coord_transform | ||
| 223 | procedure(user_function_3d_interface), optional :: userfun | ||
| 224 | procedure(user_function_3d_mat_interface), optional :: userfun_mat | ||
| 225 | integer, intent(in), optional :: n_quad_extra | ||
| 226 | |||
| 227 | integer :: s1, s2, s2_memory | ||
| 228 | |||
| 229 | 458 | type(TensorProdMat), allocatable :: M_TP_block(:, :) | |
| 230 | logical :: is_diagonal | ||
| 231 | |||
| 232 | TIMER_DECL(t_total) | ||
| 233 | TIMER_DECL(t_get_petsc) | ||
| 234 | TIMER_INIT_START(t_total, "MassMatrix::init") | ||
| 235 | |||
| 236 |
5/6✗ Branch 2 → 3 not taken.
✓ Branch 2 → 4 taken 458 times.
✓ Branch 5 → 6 taken 1374 times.
✓ Branch 5 → 10 taken 458 times.
✓ Branch 7 → 8 taken 4122 times.
✓ Branch 7 → 9 taken 1374 times.
|
5954 | allocate (M_TP_block(3, 3)) |
| 237 | |||
| 238 |
3/4✓ Branch 10 → 11 taken 278 times.
✓ Branch 10 → 12 taken 180 times.
✗ Branch 11 → 12 not taken.
✓ Branch 11 → 13 taken 278 times.
|
458 | if (.not. present(coord_transform)) then |
| 239 | 180 | is_diagonal = .not. present(userfun_mat) | |
| 240 | else | ||
| 241 |
4/4✓ Branch 13 → 14 taken 260 times.
✓ Branch 13 → 15 taken 18 times.
✓ Branch 14 → 15 taken 8 times.
✓ Branch 14 → 16 taken 252 times.
|
278 | is_diagonal = coord_transform%is_orthogonal .and. .not. present(userfun_mat) |
| 242 | end if | ||
| 243 | ! Vector valued (i.e. resulting in a block matrix) | ||
| 244 |
2/2✓ Branch 17 → 18 taken 1374 times.
✓ Branch 17 → 32 taken 458 times.
|
1832 | do s1 = 1, 3 |
| 245 |
2/2✓ Branch 19 → 20 taken 4122 times.
✓ Branch 19 → 30 taken 1374 times.
|
5954 | do s2 = 1, 3 |
| 246 |
2/2✓ Branch 20 → 21 taken 3888 times.
✓ Branch 20 → 22 taken 234 times.
|
4122 | if (is_diagonal) then |
| 247 |
2/2✓ Branch 21 → 22 taken 1296 times.
✓ Branch 21 → 29 taken 2592 times.
|
3888 | if (s1 /= s2) then |
| 248 | ! If the coordinate transformation is orthogonal, the mass matrix is | ||
| 249 | ! block diagonal | ||
| 250 | cycle | ||
| 251 | end if | ||
| 252 | s2_memory = 1 | ||
| 253 | else | ||
| 254 | s2_memory = s2 | ||
| 255 | end if | ||
| 256 | |||
| 257 | ! TODO the matrix is block symmetric, hence we need not compute all of the blocks | ||
| 258 | |||
| 259 |
4/6✓ Branch 22 → 23 taken 990 times.
✓ Branch 22 → 24 taken 540 times.
✗ Branch 23 → 24 not taken.
✓ Branch 23 → 25 taken 990 times.
✗ Branch 24 → 25 not taken.
✓ Branch 24 → 28 taken 540 times.
|
2904 | if (present(coord_transform) .or. present(userfun) .or. present(userfun_mat)) then |
| 260 |
2/2✓ Branch 25 → 26 taken 792 times.
✓ Branch 25 → 27 taken 198 times.
|
990 | if (space%m == 1) then |
| 261 | call mass_matrix(M_TP_block(s1, s2_memory), space%tp_spaces(s1), space%tp_spaces(s2), coordfun=coordfun_oneform, & | ||
| 262 | 792 | & n_quad_extra=n_quad_extra) | |
| 263 | else ! if (space%m == 2) then | ||
| 264 | call mass_matrix(M_TP_block(s1, s2_memory), space%tp_spaces(s1), space%tp_spaces(s2), coordfun=coordfun_twoform, & | ||
| 265 | 198 | & n_quad_extra=n_quad_extra) | |
| 266 | end if | ||
| 267 | else | ||
| 268 | 540 | call mass_matrix(M_TP_block(s1, s2_memory), space%tp_spaces(s1), space%tp_spaces(s2), n_quad_extra=n_quad_extra) | |
| 269 | end if | ||
| 270 | end do | ||
| 271 | end do | ||
| 272 | |||
| 273 | TIMER_INIT_START(t_get_petsc, "MassMatrix::get_petsc") | ||
| 274 |
2/2✓ Branch 33 → 34 taken 432 times.
✓ Branch 33 → 36 taken 26 times.
|
458 | if (is_diagonal) then |
| 275 | 432 | call get_petsc(this%mat, M_TP_block(:, 1)) | |
| 276 | else | ||
| 277 | 26 | call get_petsc(this%mat, M_TP_block) | |
| 278 | end if | ||
| 279 | TIMER_STOP(t_get_petsc) | ||
| 280 | |||
| 281 |
2/2✓ Branch 38 → 39 taken 1374 times.
✓ Branch 38 → 45 taken 458 times.
|
1832 | do s1 = 1, 3 |
| 282 |
2/2✓ Branch 40 → 41 taken 4122 times.
✓ Branch 40 → 43 taken 1374 times.
|
5954 | do s2 = 1, 3 |
| 283 | 5496 | call M_TP_block(s1, s2)%destroy() | |
| 284 | end do | ||
| 285 | end do | ||
| 286 | |||
| 287 |
8/12✓ Branch 46 → 47 taken 458 times.
✗ Branch 46 → 55 not taken.
✓ Branch 47 → 48 taken 4122 times.
✓ Branch 47 → 55 taken 458 times.
✗ Branch 48 → 49 not taken.
✓ Branch 48 → 50 taken 4122 times.
✗ Branch 50 → 51 not taken.
✓ Branch 50 → 52 taken 4122 times.
✓ Branch 52 → 53 taken 1530 times.
✓ Branch 52 → 54 taken 2592 times.
✗ Branch 55 → 56 not taken.
✓ Branch 55 → 57 taken 458 times.
|
4580 | deallocate (M_TP_block) |
| 288 | |||
| 289 | TIMER_STOP(t_total) | ||
| 290 | contains | ||
| 291 | |||
| 292 | 156543208 | pure real(wp) function coordfun_oneform(xp, yp, zp) | |
| 293 | use m_common, only: wp | ||
| 294 | implicit none | ||
| 295 | real(wp), intent(in) :: xp, yp, zp | ||
| 296 | |||
| 297 | integer :: m, n | ||
| 298 | real(wp) :: tmp | ||
| 299 | |||
| 300 |
2/4✓ Branch 2 → 3 taken 156543208 times.
✗ Branch 2 → 18 not taken.
✓ Branch 3 → 4 taken 156543208 times.
✗ Branch 3 → 18 not taken.
|
156543208 | if (present(coord_transform)) then |
| 301 |
2/2✓ Branch 4 → 5 taken 118127128 times.
✓ Branch 4 → 6 taken 38416080 times.
|
156543208 | if (.not. present(userfun_mat)) then |
| 302 | 118127128 | coordfun_oneform = G_matrix_inv(coord_transform, s1, s2, xp, yp, zp) | |
| 303 | else | ||
| 304 | coordfun_oneform = 0._wp | ||
| 305 |
2/2✓ Branch 7 → 8 taken 115248240 times.
✓ Branch 7 → 16 taken 38416080 times.
|
153664320 | do n = 1, 3 |
| 306 |
2/2✓ Branch 9 → 10 taken 345744720 times.
✓ Branch 9 → 14 taken 115248240 times.
|
499409040 | do m = 1, 3 |
| 307 | 345744720 | tmp = userfun_mat(m, n, xp, yp, zp) | |
| 308 |
2/2✓ Branch 11 → 12 taken 192080400 times.
✓ Branch 11 → 13 taken 153664320 times.
|
460992960 | if (tmp /= 0._wp) then |
| 309 | tmp = tmp * jacobian_matrix_inv(coord_transform, s1, m, xp, yp, zp) * & | ||
| 310 | 192080400 | jacobian_matrix_inv(coord_transform, s2, n, xp, yp, zp) | |
| 311 | |||
| 312 | 192080400 | coordfun_oneform = coordfun_oneform + tmp | |
| 313 | end if | ||
| 314 | end do | ||
| 315 | end do | ||
| 316 | end if | ||
| 317 | 156543208 | coordfun_oneform = coordfun_oneform * jacobian(coord_transform, xp, yp, zp) | |
| 318 | else | ||
| 319 | ✗ | if (.not. present(userfun_mat)) then | |
| 320 | coordfun_oneform = 1._wp | ||
| 321 | else | ||
| 322 | ✗ | coordfun_oneform = userfun_mat(s1, s2, xp, yp, zp) | |
| 323 | end if | ||
| 324 | end if | ||
| 325 | |||
| 326 |
1/2✗ Branch 20 → 21 not taken.
✓ Branch 20 → 23 taken 156543208 times.
|
156543208 | if (present(userfun)) then |
| 327 | ✗ | coordfun_oneform = coordfun_oneform * userfun(xp, yp, zp) | |
| 328 | end if | ||
| 329 | 156543208 | end function coordfun_oneform | |
| 330 | |||
| 331 | 44397075 | pure real(wp) function coordfun_twoform(xp, yp, zp) | |
| 332 | use m_common, only: wp | ||
| 333 | implicit none | ||
| 334 | real(wp), intent(in) :: xp, yp, zp | ||
| 335 | |||
| 336 | integer :: m, n | ||
| 337 | real(wp) :: tmp | ||
| 338 | |||
| 339 |
2/4✓ Branch 2 → 3 taken 44397075 times.
✗ Branch 2 → 18 not taken.
✓ Branch 3 → 4 taken 44397075 times.
✗ Branch 3 → 18 not taken.
|
44397075 | if (present(coord_transform)) then |
| 340 |
1/2✓ Branch 4 → 5 taken 44397075 times.
✗ Branch 4 → 6 not taken.
|
44397075 | if (.not. present(userfun_mat)) then |
| 341 | 44397075 | coordfun_twoform = G_matrix(coord_transform, s1, s2, xp, yp, zp) | |
| 342 | else | ||
| 343 | coordfun_twoform = 0._wp | ||
| 344 | ✗ | do n = 1, 3 | |
| 345 | ✗ | do m = 1, 3 | |
| 346 | ✗ | tmp = userfun_mat(m, n, xp, yp, zp) | |
| 347 | ✗ | if (tmp /= 0._wp) then | |
| 348 | tmp = tmp * jacobian_matrix(coord_transform, m, s1, xp, yp, zp) * & | ||
| 349 | ✗ | jacobian_matrix(coord_transform, n, s2, xp, yp, zp) | |
| 350 | |||
| 351 | ✗ | coordfun_twoform = coordfun_twoform + tmp | |
| 352 | end if | ||
| 353 | end do | ||
| 354 | end do | ||
| 355 | end if | ||
| 356 | 44397075 | coordfun_twoform = coordfun_twoform / jacobian(coord_transform, xp, yp, zp) | |
| 357 | else | ||
| 358 | ✗ | if (.not. present(userfun_mat)) then | |
| 359 | coordfun_twoform = 1._wp | ||
| 360 | else | ||
| 361 | ✗ | coordfun_twoform = userfun_mat(s1, s2, xp, yp, zp) | |
| 362 | end if | ||
| 363 | end if | ||
| 364 | |||
| 365 |
1/2✗ Branch 20 → 21 not taken.
✓ Branch 20 → 23 taken 44397075 times.
|
44397075 | if (present(userfun)) then |
| 366 | ✗ | coordfun_twoform = coordfun_twoform * userfun(xp, yp, zp) | |
| 367 | end if | ||
| 368 | 44397075 | end function coordfun_twoform | |
| 369 | end subroutine | ||
| 370 | |||
| 371 | !> @brief Initialize the matrix representation of the weak form of a differential operator for the given m-form space (gradient, | ||
| 372 | !> curl, divergence) | ||
| 373 | !> | ||
| 374 | !> @param[inout] this The DiffMatrix object to initialize | ||
| 375 | !> @param[in] space The MFormSpace for which to initialize the differential matrix | ||
| 376 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 377 | !> @param[in] mass _(optional)_ The MassMatrix corresponding to the next m-form space; if not provided, it is constructed | ||
| 378 | !> internally | ||
| 379 | ✗ | subroutine init_diff_matrix(this, space, coord_transform, mass) | |
| 380 | use m_common, only: wp | ||
| 381 | use m_mform_basis, only: next | ||
| 382 | use m_coord_transform, only: CoordTransformAbstract | ||
| 383 | implicit none | ||
| 384 | |||
| 385 | class(DiffMatrix), intent(out) :: this | ||
| 386 | type(MFormSpace), intent(in) :: space | ||
| 387 | class(CoordTransformAbstract), optional, intent(in) :: coord_transform | ||
| 388 | type(MassMatrix), optional, intent(in) :: mass | ||
| 389 | |||
| 390 | real(wp) :: EXPECTED_FILL_RATIO = 1.0_wp | ||
| 391 | integer :: ierr | ||
| 392 | |||
| 393 | Mat :: grad | ||
| 394 | ✗ | type(MassMatrix) :: mass_ | |
| 395 | |||
| 396 | TIMER_DECL(t_total) | ||
| 397 | TIMER_INIT_START(t_total, "DiffMatrix::init") | ||
| 398 | |||
| 399 | ✗ | call this%destroy() | |
| 400 | |||
| 401 | ✗ | this%space_col = space | |
| 402 | ✗ | this%space_row = next(space) | |
| 403 | |||
| 404 | ✗ | call diff_matrix_derham(grad, space, this%space_row) | |
| 405 | ✗ | if (present(mass)) then | |
| 406 | ✗ | if (mass%space_col%m /= this%space_row%m) then | |
| 407 | ✗ | error stop "DiffMatrix::init_diff_matrix: optional argument mass must be the mass matrix corresponding to next(space)" | |
| 408 | end if | ||
| 409 | ✗ | PetscCall(MatMatMult(mass%mat, grad, MAT_INITIAL_MATRIX, EXPECTED_FILL_RATIO, this%mat, ierr)) | |
| 410 | else | ||
| 411 | ✗ | call mass_%init(this%space_row, coord_transform=coord_transform) | |
| 412 | ✗ | PetscCall(MatMatMult(mass_%mat, grad, MAT_INITIAL_MATRIX, EXPECTED_FILL_RATIO, this%mat, ierr)) | |
| 413 | ✗ | call mass_%destroy() | |
| 414 | end if | ||
| 415 | ✗ | PetscCall(MatDestroy(grad, ierr)) | |
| 416 | |||
| 417 | ✗ | this%is_symmetric = .false. | |
| 418 | ✗ | this%is_spd = .false. | |
| 419 | |||
| 420 | TIMER_STOP(t_total) | ||
| 421 | ✗ | end subroutine | |
| 422 | |||
| 423 | !> @brief Initialize the matrix representation of the weak form of a differential operator for the given m-form space (divgrad, | ||
| 424 | !> curlcurl, graddiv) | ||
| 425 | !> | ||
| 426 | !> @param[inout] this The DiffDiffMatrix object to initialize | ||
| 427 | !> @param[in] space The MFormSpace for which to initialize the differential matrix | ||
| 428 | !> @param[in] coord_transform _(optional)_ The coordinate transformation associated with the m-form space | ||
| 429 | !> @param[in] mass _(optional)_ The MassMatrix corresponding to the next m-form space; if not provided, it is constructed | ||
| 430 | !> internally | ||
| 431 | 265 | subroutine init_diffdiff_matrix(this, space, coord_transform, mass) | |
| 432 | use m_common, only: wp | ||
| 433 | use m_mform_basis, only: next | ||
| 434 | use m_coord_transform, only: CoordTransformAbstract | ||
| 435 | implicit none | ||
| 436 | |||
| 437 | class(DiffDiffMatrix), intent(out) :: this | ||
| 438 | type(MFormSpace), intent(in) :: space | ||
| 439 | class(CoordTransformAbstract), optional, intent(in) :: coord_transform | ||
| 440 | type(MassMatrix), optional, intent(in) :: mass | ||
| 441 | |||
| 442 | real(wp) :: EXPECTED_FILL_RATIO = 1.0_wp | ||
| 443 | integer :: ierr | ||
| 444 | |||
| 445 | Mat :: grad | ||
| 446 | 265 | type(MassMatrix) :: mass_ | |
| 447 | 265 | type(MFormSpace) :: next_space | |
| 448 | |||
| 449 | TIMER_DECL(t_total) | ||
| 450 | TIMER_INIT_START(t_total, "DiffDiffMatrix::init") | ||
| 451 | |||
| 452 | 265 | call this%destroy() | |
| 453 | |||
| 454 |
4/8✓ Branch 6 → 7 taken 265 times.
✗ Branch 6 → 10 not taken.
✓ Branch 7 → 8 taken 265 times.
✗ Branch 7 → 9 not taken.
✓ Branch 10 → 11 taken 265 times.
✗ Branch 10 → 13 not taken.
✗ Branch 11 → 12 not taken.
✓ Branch 11 → 13 taken 265 times.
|
265 | this%space_col = space |
| 455 |
4/8✓ Branch 13 → 14 taken 265 times.
✗ Branch 13 → 17 not taken.
✓ Branch 14 → 15 taken 265 times.
✗ Branch 14 → 16 not taken.
✓ Branch 17 → 18 taken 265 times.
✗ Branch 17 → 20 not taken.
✗ Branch 18 → 19 not taken.
✓ Branch 18 → 20 taken 265 times.
|
265 | this%space_row = space |
| 456 | |||
| 457 | 265 | next_space = next(space) | |
| 458 | |||
| 459 | 265 | call diff_matrix_derham(grad, space, next_space) | |
| 460 |
2/2✓ Branch 22 → 23 taken 38 times.
✓ Branch 22 → 28 taken 227 times.
|
265 | if (present(mass)) then |
| 461 |
1/2✗ Branch 23 → 24 not taken.
✓ Branch 23 → 25 taken 38 times.
|
38 | if (mass%space_col%m /= next_space%m) then |
| 462 | error stop "DiffDiffMatrix::init_diffdiff_matrix: optional argument mass must be the mass matrix corresponding to"// & | ||
| 463 | ✗ | &" next(space)" | |
| 464 | end if | ||
| 465 |
1/2✗ Branch 26 → 27 not taken.
✓ Branch 26 → 34 taken 38 times.
|
38 | PetscCall(MatPtAP(mass%mat, grad, MAT_INITIAL_MATRIX, EXPECTED_FILL_RATIO, this%mat, ierr)) |
| 466 | else | ||
| 467 | 227 | call mass_%init(next_space, coord_transform=coord_transform) | |
| 468 |
1/2✗ Branch 30 → 31 not taken.
✓ Branch 30 → 32 taken 227 times.
|
227 | PetscCall(MatPtAP(mass_%mat, grad, MAT_INITIAL_MATRIX, EXPECTED_FILL_RATIO, this%mat, ierr)) |
| 469 | 227 | call mass_%destroy() | |
| 470 | end if | ||
| 471 |
1/2✗ Branch 35 → 36 not taken.
✓ Branch 35 → 37 taken 265 times.
|
265 | PetscCall(MatDestroy(grad, ierr)) |
| 472 | |||
| 473 | 265 | this%is_symmetric = .true. | |
| 474 | 265 | this%is_spd = .true. | |
| 475 | |||
| 476 | TIMER_STOP(t_total) | ||
| 477 |
6/14✓ Branch 2 → 3 taken 265 times.
✗ Branch 2 → 5 not taken.
✓ Branch 37 → 38 taken 265 times.
✗ Branch 37 → 39 not taken.
✓ Branch 39 → 40 taken 227 times.
✓ Branch 39 → 41 taken 38 times.
✓ Branch 41 → 42 taken 227 times.
✓ Branch 41 → 54 taken 38 times.
✗ Branch 43 → 44 not taken.
✗ Branch 43 → 45 not taken.
✗ Branch 47 → 48 not taken.
✗ Branch 47 → 49 not taken.
✗ Branch 51 → 52 not taken.
✗ Branch 51 → 55 not taken.
|
530 | end subroutine |
| 478 | |||
| 479 | !> @brief Initialize a zero matrix for the given m-form spaces | ||
| 480 | !> | ||
| 481 | !> @param[inout] this The ZeroMatrix object to initialize | ||
| 482 | !> @param[in] space_row The m-form space corresponding to the rows of the matrix | ||
| 483 | !> @param[in] space_col _(optional)_ The m-form space corresponding to the columns of the matrix | ||
| 484 | !> | ||
| 485 | !> @note If space_row is not provided, a square zero matrix on space_col is constructed | ||
| 486 | 8 | subroutine init_zero_matrix(this, space_row, space_col) | |
| 487 | use m_mform_basis, only: size | ||
| 488 | implicit none | ||
| 489 | |||
| 490 | class(ZeroMatrix), intent(inout) :: this | ||
| 491 | type(MFormSpace), intent(in) :: space_row | ||
| 492 | type(MFormSpace), optional, intent(in) :: space_col | ||
| 493 | |||
| 494 | integer :: ierr | ||
| 495 | integer :: local_nz_per_row, nonlocal_nz_per_row | ||
| 496 | |||
| 497 | TIMER_DECL(t_total) | ||
| 498 | TIMER_INIT_START(t_total, "ZeroMatrix::init") | ||
| 499 | |||
| 500 | 8 | call this%destroy() | |
| 501 | |||
| 502 |
4/8✓ Branch 3 → 4 taken 8 times.
✗ Branch 3 → 7 not taken.
✓ Branch 4 → 5 taken 8 times.
✗ Branch 4 → 6 not taken.
✓ Branch 7 → 8 taken 8 times.
✗ Branch 7 → 10 not taken.
✗ Branch 8 → 9 not taken.
✓ Branch 8 → 10 taken 8 times.
|
8 | this%space_row = space_row |
| 503 |
1/2✓ Branch 10 → 11 taken 8 times.
✗ Branch 10 → 19 not taken.
|
8 | if (present(space_col)) then |
| 504 |
4/8✓ Branch 11 → 12 taken 8 times.
✗ Branch 11 → 15 not taken.
✓ Branch 12 → 13 taken 8 times.
✗ Branch 12 → 14 not taken.
✓ Branch 15 → 16 taken 8 times.
✗ Branch 15 → 18 not taken.
✗ Branch 16 → 17 not taken.
✓ Branch 16 → 18 taken 8 times.
|
8 | this%space_col = space_col |
| 505 | else | ||
| 506 | ✗ | this%space_col = space_row | |
| 507 | end if | ||
| 508 | |||
| 509 | 8 | local_nz_per_row = 0 | |
| 510 | 8 | nonlocal_nz_per_row = 0 | |
| 511 | call MatCreateAIJ(space_col%tp_spaces(1)%domain%comm, size(this%space_row, my_rank=.true.), & | ||
| 512 | & size(this%space_col, my_rank=.true.), size(this%space_row), size(this%space_col), & | ||
| 513 | 8 | & local_nz_per_row, PETSC_NULL_INTEGER_ARRAY, nonlocal_nz_per_row, PETSC_NULL_INTEGER_ARRAY, this%mat, ierr) | |
| 514 |
1/2✗ Branch 28 → 29 not taken.
✓ Branch 28 → 30 taken 8 times.
|
8 | CHKERRQ(ierr) |
| 515 | |||
| 516 |
1/2✗ Branch 31 → 32 not taken.
✓ Branch 31 → 33 taken 8 times.
|
8 | PetscCall(MatSetFromOptions(this%mat, ierr)) |
| 517 |
1/2✗ Branch 34 → 35 not taken.
✓ Branch 34 → 36 taken 8 times.
|
8 | PetscCall(MatSetUp(this%mat, ierr)) |
| 518 |
1/2✗ Branch 37 → 38 not taken.
✓ Branch 37 → 39 taken 8 times.
|
8 | PetscCall(MatAssemblyBegin(this%mat, MAT_FINAL_ASSEMBLY, ierr)) |
| 519 |
1/2✗ Branch 40 → 41 not taken.
✓ Branch 40 → 42 taken 8 times.
|
8 | PetscCall(MatAssemblyEnd(this%mat, MAT_FINAL_ASSEMBLY, ierr)) |
| 520 | |||
| 521 | 8 | this%is_symmetric = (size(this%space_row) == size(this%space_col)) | |
| 522 | 8 | this%is_spd = .false. | |
| 523 | |||
| 524 | TIMER_STOP(t_total) | ||
| 525 | end subroutine init_zero_matrix | ||
| 526 | |||
| 527 | !> @brief Coefficient for the imposing a linear constraint on the mass matrix by means of a projection operator | ||
| 528 | !> | ||
| 529 | !> That is, when mode=CONSTRAINT_PROJECTION is used, the following matrix is constructed: | ||
| 530 | !> P^T * M * P + \gamma * (I - P)^T M * (I - P), | ||
| 531 | !> where P is the projection operator and M is the mass matrix. | ||
| 532 | !> This function provides the coefficient \gamma. | ||
| 533 | !> | ||
| 534 | !> @param[in] this The MassMatrix object for which to compute the coefficient | ||
| 535 | !> | ||
| 536 | !> @return The coefficient \gamma | ||
| 537 | 105 | pure function constraint_coeff_mass(this) result(ans) | |
| 538 | use m_common, only: wp | ||
| 539 | |||
| 540 | class(MassMatrix), intent(in) :: this | ||
| 541 | |||
| 542 | real(wp) :: ans | ||
| 543 | |||
| 544 | ans = 1.0_wp | ||
| 545 | 105 | end function constraint_coeff_mass | |
| 546 | |||
| 547 | !> @brief Coefficient for the imposing a linear constraint on the differential matrix by means of a projection operator | ||
| 548 | !> | ||
| 549 | !> That is, when mode=CONSTRAINT_PROJECTION is used, the following matrix is constructed: | ||
| 550 | !> P^T * L * D + \gamma * (I - P)^T M * (I - P), | ||
| 551 | !> where P is the projection operator, L is the divgrad/curlcurl/graddiv matrix. | ||
| 552 | !> This function provides the coefficient \gamma. | ||
| 553 | !> | ||
| 554 | !> @param[in] this The DiffDiffMatrix object for which to compute the coefficient | ||
| 555 | !> | ||
| 556 | !> @return The coefficient \gamma | ||
| 557 | 257 | pure function constraint_coeff_diffdiff(this) result(ans) | |
| 558 | use m_common, only: wp | ||
| 559 | |||
| 560 | class(DiffDiffMatrix), intent(in) :: this | ||
| 561 | |||
| 562 | real(wp) :: ans | ||
| 563 | |||
| 564 | ! TODO is this good? at least now it's dimensionally consistent | ||
| 565 | ans = 1._wp / ((this%space_col%tp_spaces(1)%spaces(1)%h * this%space_col%tp_spaces(1)%spaces(2)%h * & | ||
| 566 | 257 | & this%space_col%tp_spaces(1)%spaces(3)%h)**(1._wp / 3._wp)) | |
| 567 | 257 | end function constraint_coeff_diffdiff | |
| 568 | |||
| 569 | !> @brief This function is meaningless for the DiffDiffMatrix, but is required by the AbstractMatrix interface | ||
| 570 | ✗ | pure function constraint_coeff_diff(this) result(ans) | |
| 571 | use m_common, only: wp | ||
| 572 | |||
| 573 | class(DiffMatrix), intent(in) :: this | ||
| 574 | |||
| 575 | real(wp) :: ans | ||
| 576 | |||
| 577 | ! NOTE: Never used but AbstractMatrix requires it | ||
| 578 | ! TODO maybe have intermediate type 'AbstractSquareMatrix'? | ||
| 579 | ans = 1._wp | ||
| 580 | ✗ | end function constraint_coeff_diff | |
| 581 | |||
| 582 | !> @brief This function is meaningless for the ZeroMatrix, but is required by the AbstractMatrix interface | ||
| 583 | 16 | pure function constraint_coeff_zero(this) result(ans) | |
| 584 | use m_common, only: wp | ||
| 585 | |||
| 586 | class(ZeroMatrix), intent(in) :: this | ||
| 587 | |||
| 588 | real(wp) :: ans | ||
| 589 | |||
| 590 | ans = 1._wp | ||
| 591 | 16 | end function constraint_coeff_zero | |
| 592 | |||
| 593 | !> @brief Destroy the AbstractMatrix | ||
| 594 | !> | ||
| 595 | !> @param[inout] this The AbstractMatrix object to destroy | ||
| 596 | 1534 | subroutine destroy_abstract_matrix(this) | |
| 597 | implicit none | ||
| 598 | |||
| 599 | class(AbstractMatrix), intent(inout) :: this | ||
| 600 | integer :: ierr | ||
| 601 | |||
| 602 |
1/2✗ Branch 3 → 4 not taken.
✓ Branch 3 → 6 taken 1534 times.
|
1534 | PetscCall(MatDestroy(this%mat, ierr)) |
| 603 | end subroutine destroy_abstract_matrix | ||
| 604 | |||
| 605 | !> @brief Create the matrix representation of the weak form of a differential operator grad/curl/div in the DeRham sequence | ||
| 606 | !> | ||
| 607 | !> @param[out] diff The matrix to create | ||
| 608 | !> @param[in] space The MFormSpace for which to create the differential matrix | ||
| 609 | !> @param[in] next_space The MFormSpace for the next m-form space in the DeRham sequence (i.e., next_space=next(space)) | ||
| 610 | 265 | subroutine diff_matrix_derham(diff, space, next_space) | |
| 611 | implicit none | ||
| 612 | |||
| 613 | Mat, intent(out) :: diff | ||
| 614 | type(MFormSpace), intent(in) :: space, next_space | ||
| 615 | |||
| 616 |
1/2✗ Branch 2 → 3 not taken.
✓ Branch 2 → 4 taken 265 times.
|
265 | if (space%m + 1 /= next_space%m) then |
| 617 | ✗ | error stop "diff_matrix_derham: space and next_space must be consecutive in the DeRham sequence" | |
| 618 | end if | ||
| 619 | |||
| 620 | 237 | select case (space%m) | |
| 621 | case (0) | ||
| 622 | 265 | call gradient_matrix(diff, space, next_space) | |
| 623 | case (1) | ||
| 624 | 28 | call curl_matrix(diff, space, next_space) | |
| 625 | case (2) | ||
| 626 |
2/4✓ Branch 4 → 5 taken 237 times.
✓ Branch 4 → 7 taken 28 times.
✗ Branch 4 → 9 not taken.
✗ Branch 4 → 10 not taken.
|
265 | call divergence_matrix(diff, space, next_space) |
| 627 | end select | ||
| 628 | 265 | end subroutine diff_matrix_derham | |
| 629 | |||
| 630 | !> @brief Create the matrix representation of the weak form of the gradient operator in the DeRham sequence | ||
| 631 | !> | ||
| 632 | !> @param[out] grad The matrix to create | ||
| 633 | !> @param[in] space The MFormSpace for which to create the gradient matrix | ||
| 634 | !> @param[in] next_space The MFormSpace for the next m-form space in the DeRham sequence (i.e., next_space=next(space)) | ||
| 635 | 262 | subroutine gradient_matrix(grad, space, next_space) | |
| 636 | use m_tensorprod, only: size | ||
| 637 | use m_mform_basis, only: size | ||
| 638 | implicit none | ||
| 639 | |||
| 640 | Mat, intent(out) :: grad | ||
| 641 | type(MFormSpace), intent(in) :: space, next_space | ||
| 642 | |||
| 643 | integer :: ierr, local_nz_per_row, nonlocal_nz_per_row | ||
| 644 | integer :: blk | ||
| 645 | |||
| 646 | TIMER_DECL(t_total) | ||
| 647 | TIMER_INIT_START(t_total, "MFormSpace::gradient_matrix") | ||
| 648 | |||
| 649 | 262 | local_nz_per_row = 2 | |
| 650 | 262 | nonlocal_nz_per_row = 2 | |
| 651 | call MatCreateAIJ(space%tp_spaces(1)%domain%comm, size(next_space, my_rank=.true.), size(space, my_rank=.true.), & | ||
| 652 | & size(next_space), size(space), & | ||
| 653 | & local_nz_per_row, PETSC_NULL_INTEGER_ARRAY, nonlocal_nz_per_row, & | ||
| 654 |
1/2✗ Branch 3 → 4 not taken.
✓ Branch 3 → 5 taken 262 times.
|
262 | & PETSC_NULL_INTEGER_ARRAY, grad, ierr); CHKERRQ(ierr) |
| 655 | |||
| 656 |
1/2✗ Branch 6 → 7 not taken.
✓ Branch 6 → 8 taken 262 times.
|
262 | PetscCall(MatSetFromOptions(grad, ierr)) |
| 657 |
1/2✗ Branch 9 → 10 not taken.
✓ Branch 9 → 11 taken 262 times.
|
262 | PetscCall(MatSetUp(grad, ierr)) |
| 658 | |||
| 659 |
2/2✓ Branch 12 → 13 taken 554 times.
✓ Branch 12 → 15 taken 262 times.
|
816 | do blk = 1, next_space%dimensionality |
| 660 | 816 | call insert_diff_matrix(grad, next_space%tp_spaces(blk), space%tp_spaces(1), blk) | |
| 661 | end do | ||
| 662 | |||
| 663 |
1/2✗ Branch 17 → 18 not taken.
✓ Branch 17 → 19 taken 262 times.
|
262 | PetscCall(MatAssemblyBegin(grad, MAT_FINAL_ASSEMBLY, ierr)) |
| 664 |
1/2✗ Branch 20 → 21 not taken.
✓ Branch 20 → 22 taken 262 times.
|
262 | PetscCall(MatAssemblyEnd(grad, MAT_FINAL_ASSEMBLY, ierr)) |
| 665 | |||
| 666 | TIMER_STOP(t_total) | ||
| 667 | end subroutine gradient_matrix | ||
| 668 | |||
| 669 | !> @brief Create the matrix representation of the weak form of the curl operator in the DeRham sequence | ||
| 670 | !> | ||
| 671 | !> @param[out] curl The matrix to create | ||
| 672 | !> @param[in] space The MFormSpace for which to create the curl matrix | ||
| 673 | !> @param[in] next_space The MFormSpace for the next m-form space in the DeRham sequence (i.e., next_space=next(space)) | ||
| 674 | 30 | subroutine curl_matrix(curl, space, next_space) | |
| 675 | use m_tensorprod, only: size | ||
| 676 | use m_common, only: levi_civita, wp | ||
| 677 | use m_mform_basis, only: size | ||
| 678 | implicit none | ||
| 679 | |||
| 680 | Mat, intent(out) :: curl | ||
| 681 | type(MFormSpace), intent(in) :: space, next_space | ||
| 682 | |||
| 683 | integer :: ierr, local_nz_per_row, nonlocal_nz_per_row | ||
| 684 | integer :: row_blk, col_blk, dir, lc | ||
| 685 | |||
| 686 | TIMER_DECL(t_total) | ||
| 687 | TIMER_INIT_START(t_total, "MFormSpace::curl_matrix") | ||
| 688 | |||
| 689 | 30 | local_nz_per_row = 4 | |
| 690 | 30 | nonlocal_nz_per_row = 4 | |
| 691 | call MatCreateAIJ(space%tp_spaces(1)%domain%comm, size(next_space, my_rank=.true.), size(space, my_rank=.true.), & | ||
| 692 | & size(next_space), size(space), & | ||
| 693 | & local_nz_per_row, PETSC_NULL_INTEGER_ARRAY, nonlocal_nz_per_row, & | ||
| 694 |
1/2✗ Branch 3 → 4 not taken.
✓ Branch 3 → 5 taken 30 times.
|
30 | & PETSC_NULL_INTEGER_ARRAY, curl, ierr); CHKERRQ(ierr) |
| 695 | |||
| 696 |
1/2✗ Branch 6 → 7 not taken.
✓ Branch 6 → 8 taken 30 times.
|
30 | PetscCall(MatSetFromOptions(curl, ierr)) |
| 697 |
1/2✗ Branch 9 → 10 not taken.
✓ Branch 9 → 11 taken 30 times.
|
30 | PetscCall(MatSetUp(curl, ierr)) |
| 698 | |||
| 699 |
2/2✓ Branch 12 → 13 taken 90 times.
✓ Branch 12 → 25 taken 30 times.
|
120 | do row_blk = 1, 3 |
| 700 |
2/2✓ Branch 14 → 15 taken 270 times.
✓ Branch 14 → 23 taken 90 times.
|
390 | do col_blk = 1, 3 |
| 701 |
2/2✓ Branch 16 → 17 taken 747 times.
✓ Branch 16 → 21 taken 270 times.
|
1107 | do dir = 1, space%dimensionality |
| 702 | 747 | lc = levi_civita(row_blk, dir, col_blk) | |
| 703 |
2/2✓ Branch 17 → 18 taken 166 times.
✓ Branch 17 → 20 taken 581 times.
|
1017 | if (lc /= 0) then |
| 704 | call insert_diff_matrix(curl, next_space%tp_spaces(row_blk), space%tp_spaces(col_blk), & | ||
| 705 | 166 | dir, multiply_by=real(lc, kind=wp)) | |
| 706 | end if | ||
| 707 | end do | ||
| 708 | end do | ||
| 709 | end do | ||
| 710 | |||
| 711 |
1/2✗ Branch 27 → 28 not taken.
✓ Branch 27 → 29 taken 30 times.
|
30 | PetscCall(MatAssemblyBegin(curl, MAT_FINAL_ASSEMBLY, ierr)) |
| 712 |
1/2✗ Branch 30 → 31 not taken.
✓ Branch 30 → 32 taken 30 times.
|
30 | PetscCall(MatAssemblyEnd(curl, MAT_FINAL_ASSEMBLY, ierr)) |
| 713 | |||
| 714 | TIMER_STOP(t_total) | ||
| 715 | end subroutine curl_matrix | ||
| 716 | |||
| 717 | !> @brief Create the matrix representation of the weak form of the divergence operator in the DeRham sequence | ||
| 718 | !> | ||
| 719 | !> @param[out] div The matrix to create | ||
| 720 | !> @param[in] space The MFormSpace for which to create the divergence matrix | ||
| 721 | !> @param[in] next_space The MFormSpace for the next m-form space in the DeRham sequence (i.e., next_space=next(space)) | ||
| 722 | ✗ | subroutine divergence_matrix(div, space, next_space) | |
| 723 | use m_tensorprod, only: size | ||
| 724 | use m_mform_basis, only: size | ||
| 725 | implicit none | ||
| 726 | |||
| 727 | Mat, intent(out) :: div | ||
| 728 | type(MFormSpace), intent(in) :: space, next_space | ||
| 729 | |||
| 730 | integer :: ierr, local_nz_per_row, nonlocal_nz_per_row | ||
| 731 | integer :: blk | ||
| 732 | |||
| 733 | TIMER_DECL(t_total) | ||
| 734 | TIMER_INIT_START(t_total, "MFormSpace::divergence_matrix") | ||
| 735 | |||
| 736 | ✗ | local_nz_per_row = 6 | |
| 737 | ✗ | nonlocal_nz_per_row = 6 | |
| 738 | call MatCreateAIJ(space%tp_spaces(1)%domain%comm, size(next_space, my_rank=.true.), size(space, my_rank=.true.), & | ||
| 739 | & size(next_space), size(space), & | ||
| 740 | & local_nz_per_row, PETSC_NULL_INTEGER_ARRAY, nonlocal_nz_per_row, & | ||
| 741 | ✗ | & PETSC_NULL_INTEGER_ARRAY, div, ierr); CHKERRQ(ierr) | |
| 742 | |||
| 743 | ✗ | PetscCall(MatSetFromOptions(div, ierr)) | |
| 744 | ✗ | PetscCall(MatSetUp(div, ierr)) | |
| 745 | |||
| 746 | ✗ | do blk = 1, space%dimensionality | |
| 747 | ✗ | call insert_diff_matrix(div, next_space%tp_spaces(1), space%tp_spaces(blk), blk) | |
| 748 | end do | ||
| 749 | |||
| 750 | ✗ | PetscCall(MatAssemblyBegin(div, MAT_FINAL_ASSEMBLY, ierr)) | |
| 751 | ✗ | PetscCall(MatAssemblyEnd(div, MAT_FINAL_ASSEMBLY, ierr)) | |
| 752 | |||
| 753 | TIMER_STOP(t_total) | ||
| 754 | end subroutine divergence_matrix | ||
| 755 | |||
| 756 | !> @brief Insert the finite difference matrix for the derivative in the given direction into the matrix | ||
| 757 | !> | ||
| 758 | !> @param[inout] mat The matrix to insert the finite difference matrix into | ||
| 759 | !> @param[in] row_space The MFormSpace corresponding to the rows of the matrix | ||
| 760 | !> @param[in] col_space The MFormSpace corresponding to the columns of the matrix | ||
| 761 | !> @param[in] diff_dir The direction of the derivative (1 for 'x', 2 for 'y', 3 for 'z') | ||
| 762 | !> @param[in] multiply_by _(optional)_ A scaling factor to apply to the finite difference matrix | ||
| 763 | !> | ||
| 764 | !> @note See also the 'add_difference' subroutine in m_mform_basis.F90 | ||
| 765 | 720 | subroutine insert_diff_matrix(mat, row_space, col_space, diff_dir, multiply_by) | |
| 766 | use m_common, only: wp | ||
| 767 | use m_tensorprod, only: size, TensorProdSpace, cartesian_to_linear | ||
| 768 | |||
| 769 | implicit none | ||
| 770 | |||
| 771 | Mat, intent(inout) :: mat | ||
| 772 | type(TensorProdSpace), intent(in) :: row_space, col_space ! NOTE: row_space = derivative(col_space, diff_dir) | ||
| 773 | integer, intent(in) :: diff_dir | ||
| 774 | real(wp), optional, intent(in) :: multiply_by | ||
| 775 | |||
| 776 | integer :: ierr | ||
| 777 | integer :: ir, jr, kr, ic, jc, kc | ||
| 778 | integer :: row, col_c, col_r | ||
| 779 | |||
| 780 | integer :: col_sze1, col_sze2, col_sze3 | ||
| 781 | |||
| 782 | real(wp) :: multiply_by_ | ||
| 783 | |||
| 784 |
2/2✓ Branch 2 → 3 taken 166 times.
✓ Branch 2 → 4 taken 554 times.
|
720 | if (present(multiply_by)) then |
| 785 | 166 | multiply_by_ = multiply_by | |
| 786 | else | ||
| 787 | multiply_by_ = 1._wp | ||
| 788 | end if | ||
| 789 | |||
| 790 | 720 | col_sze1 = size(col_space, 1) | |
| 791 | 720 | col_sze2 = size(col_space, 2) | |
| 792 | 720 | col_sze3 = size(col_space, 3) | |
| 793 | |||
| 794 |
2/2✓ Branch 5 → 6 taken 5233 times.
✓ Branch 5 → 35 taken 720 times.
|
5953 | do kc = row_space%rank_resp_bspline%k0, row_space%rank_resp_bspline%k1 |
| 795 |
2/2✓ Branch 7 → 8 taken 69643 times.
✓ Branch 7 → 33 taken 5233 times.
|
75596 | do jc = row_space%rank_resp_bspline%j0, row_space%rank_resp_bspline%j1 |
| 796 |
2/2✓ Branch 9 → 10 taken 1127380 times.
✓ Branch 9 → 31 taken 69643 times.
|
1202256 | do ic = row_space%rank_resp_bspline%i0, row_space%rank_resp_bspline%i1 |
| 797 | 1127380 | call cartesian_to_linear(row_space, row, ic, jc, kc, is_on_my_responsible_subdomain=.true.) | |
| 798 | 1127380 | call cartesian_to_linear(col_space, col_c, ic, jc, kc) | |
| 799 | |||
| 800 | 1127380 | ir = ic | |
| 801 | 1127380 | jr = jc | |
| 802 | 1127380 | kr = kc | |
| 803 | 1536026 | select case (diff_dir) | |
| 804 | case (1) | ||
| 805 | 408646 | ir = mod(ic + 1, col_sze1) | |
| 806 | case (2) | ||
| 807 | 411638 | jr = mod(jc + 1, col_sze2) | |
| 808 | case (3) | ||
| 809 |
3/4✓ Branch 12 → 13 taken 408646 times.
✓ Branch 12 → 14 taken 411638 times.
✓ Branch 12 → 15 taken 307096 times.
✗ Branch 12 → 16 not taken.
|
1127380 | kr = mod(kc + 1, col_sze3) |
| 810 | end select | ||
| 811 | |||
| 812 | 1127380 | call cartesian_to_linear(col_space, col_r, ir, jr, kr) | |
| 813 |
7/8✓ Branch 18 → 19 taken 1127380 times.
✓ Branch 18 → 20 taken 1127380 times.
✓ Branch 21 → 22 taken 2254760 times.
✓ Branch 21 → 23 taken 1127380 times.
✓ Branch 24 → 25 taken 2254760 times.
✓ Branch 24 → 26 taken 1127380 times.
✗ Branch 27 → 28 not taken.
✓ Branch 27 → 30 taken 1127380 times.
|
7961303 | PetscCall(MatSetValues(mat, 1, [row], 2, [col_c, col_r], [-multiply_by_, multiply_by_], INSERT_VALUES, ierr)) |
| 814 | end do | ||
| 815 | end do | ||
| 816 | end do | ||
| 817 | |||
| 818 | end subroutine insert_diff_matrix | ||
| 819 | |||
| 820 | !> @brief Transpose the matrix | ||
| 821 | !> | ||
| 822 | !> @param[out] mat_t The transposed matrix | ||
| 823 | !> @param[in] mat The matrix to transpose | ||
| 824 |
1/2✓ Branch 2 → 3 taken 8 times.
✗ Branch 2 → 5 not taken.
|
8 | subroutine transpose_abstract_matrix(mat_t, mat) |
| 825 | use m_common, only: wp | ||
| 826 | implicit none | ||
| 827 | |||
| 828 | class(AbstractMatrix), intent(in) :: mat | ||
| 829 | class(AbstractMatrix), intent(out) :: mat_t | ||
| 830 | |||
| 831 | integer :: ierr | ||
| 832 | |||
| 833 | TIMER_DECL(t_total) | ||
| 834 | TIMER_INIT_START(t_total, "AbstractMatrix::transpose") | ||
| 835 | |||
| 836 | 8 | call mat_t%destroy() | |
| 837 | |||
| 838 |
4/8✓ Branch 6 → 7 taken 8 times.
✗ Branch 6 → 10 not taken.
✓ Branch 7 → 8 taken 8 times.
✗ Branch 7 → 9 not taken.
✓ Branch 10 → 11 taken 8 times.
✗ Branch 10 → 13 not taken.
✗ Branch 11 → 12 not taken.
✓ Branch 11 → 13 taken 8 times.
|
8 | mat_t%space_col = mat%space_row |
| 839 |
4/8✓ Branch 13 → 14 taken 8 times.
✗ Branch 13 → 17 not taken.
✓ Branch 14 → 15 taken 8 times.
✗ Branch 14 → 16 not taken.
✓ Branch 17 → 18 taken 8 times.
✗ Branch 17 → 20 not taken.
✗ Branch 18 → 19 not taken.
✓ Branch 18 → 20 taken 8 times.
|
8 | mat_t%space_row = mat%space_col |
| 840 | |||
| 841 |
1/2✗ Branch 20 → 21 not taken.
✓ Branch 20 → 24 taken 8 times.
|
8 | if (mat%is_symmetric) then |
| 842 | ✗ | PetscCall(MatDuplicate(mat%mat, MAT_COPY_VALUES, mat_t%mat, ierr)) | |
| 843 | else | ||
| 844 |
1/2✗ Branch 25 → 26 not taken.
✓ Branch 25 → 27 taken 8 times.
|
8 | PetscCall(MatTranspose(mat%mat, MAT_INITIAL_MATRIX, mat_t%mat, ierr)) |
| 845 | end if | ||
| 846 | |||
| 847 | 8 | mat_t%is_symmetric = mat%is_symmetric | |
| 848 | 8 | mat_t%is_spd = mat%is_spd | |
| 849 | |||
| 850 | TIMER_STOP(t_total) | ||
| 851 | end subroutine transpose_abstract_matrix | ||
| 852 | |||
| 853 | !> @brief Matrix-vector product for the mass matrix | ||
| 854 | !> | ||
| 855 | !> @param[in] this The MassMatrix to use for the matrix-vector product | ||
| 856 | !> @param[in] vin The input MFormFun to multiply by the mass matrix | ||
| 857 | !> @param[out] vout The output MFormFun to store the result of the matrix-vector product | ||
| 858 | ✗ | subroutine matvec_mass(this, vout, vin) | |
| 859 | use m_common, only: wp | ||
| 860 | use m_mform_basis, only: MFormFun, get_petsc | ||
| 861 | implicit none | ||
| 862 | |||
| 863 | class(MassMatrix), intent(in) :: this | ||
| 864 | type(MFormFun), intent(in) :: vin | ||
| 865 | type(MFormFun), intent(inout) :: vout | ||
| 866 | |||
| 867 | integer :: ierr | ||
| 868 | |||
| 869 | Vec :: vec_in, vec_out | ||
| 870 | |||
| 871 | TIMER_DECL(t_total) | ||
| 872 | TIMER_INIT_START(t_total, "MassMatrix::matvec") | ||
| 873 | |||
| 874 | ! TODO make it mass matrix specific and avoid petsc | ||
| 875 | ✗ | call get_petsc(vec_in, vin) | |
| 876 | |||
| 877 | ✗ | PetscCall(VecDuplicate(vec_in, vec_out, ierr)) | |
| 878 | ✗ | PetscCall(MatMult(this%mat, vec_in, vec_out, ierr)) | |
| 879 | |||
| 880 | ✗ | call vout%reset(vec_out) | |
| 881 | |||
| 882 | TIMER_STOP(t_total) | ||
| 883 | ✗ | end subroutine matvec_mass | |
| 884 |
20/86__m_mform_matrix_MOD___copy_m_mform_matrix_Diffdiffmatrix:
✓ Branch 2 → 3 taken 16 times.
✗ Branch 2 → 9 not taken.
✓ Branch 3 → 4 taken 16 times.
✗ Branch 3 → 5 not taken.
✓ Branch 6 → 7 taken 16 times.
✗ Branch 6 → 8 not taken.
__m_mform_matrix_MOD___copy_m_mform_matrix_Diffmatrix:
✗ Branch 2 → 3 not taken.
✗ Branch 2 → 9 not taken.
✗ Branch 3 → 4 not taken.
✗ Branch 3 → 5 not taken.
✗ Branch 6 → 7 not taken.
✗ Branch 6 → 8 not taken.
__m_mform_matrix_MOD___copy_m_mform_matrix_Massmatrix:
✗ Branch 2 → 3 not taken.
✗ Branch 2 → 9 not taken.
✗ Branch 3 → 4 not taken.
✗ Branch 3 → 5 not taken.
✗ Branch 6 → 7 not taken.
✗ Branch 6 → 8 not taken.
__m_mform_matrix_MOD___copy_m_mform_matrix_Zeromatrix:
✓ Branch 2 → 3 taken 16 times.
✗ Branch 2 → 9 not taken.
✓ Branch 3 → 4 taken 16 times.
✗ Branch 3 → 5 not taken.
✓ Branch 6 → 7 taken 16 times.
✗ Branch 6 → 8 not taken.
__m_mform_matrix_MOD___final_m_mform_matrix_Abstractmatrix:
✗ Branch 3 → 4 not taken.
✗ Branch 3 → 7 not taken.
✗ Branch 4 → 5 not taken.
✗ Branch 4 → 6 not taken.
✗ Branch 8 → 9 not taken.
✗ Branch 8 → 19 not taken.
✗ Branch 10 → 11 not taken.
✗ Branch 10 → 12 not taken.
✗ Branch 12 → 13 not taken.
✗ Branch 12 → 16 not taken.
✗ Branch 13 → 14 not taken.
✗ Branch 13 → 15 not taken.
✗ Branch 16 → 17 not taken.
✗ Branch 16 → 18 not taken.
__m_mform_matrix_MOD___final_m_mform_matrix_Diffdiffmatrix:
✗ Branch 3 → 4 not taken.
✓ Branch 3 → 7 taken 281 times.
✗ Branch 4 → 5 not taken.
✗ Branch 4 → 6 not taken.
✓ Branch 8 → 9 taken 281 times.
✓ Branch 8 → 17 taken 281 times.
✗ Branch 10 → 11 not taken.
✓ Branch 10 → 12 taken 281 times.
✓ Branch 12 → 13 taken 281 times.
✗ Branch 12 → 16 not taken.
✓ Branch 13 → 14 taken 120 times.
✓ Branch 13 → 15 taken 161 times.
__m_mform_matrix_MOD___final_m_mform_matrix_Diffmatrix:
✗ Branch 3 → 4 not taken.
✗ Branch 3 → 7 not taken.
✗ Branch 4 → 5 not taken.
✗ Branch 4 → 6 not taken.
✗ Branch 8 → 9 not taken.
✗ Branch 8 → 17 not taken.
✗ Branch 10 → 11 not taken.
✗ Branch 10 → 12 not taken.
✗ Branch 12 → 13 not taken.
✗ Branch 12 → 16 not taken.
✗ Branch 13 → 14 not taken.
✗ Branch 13 → 15 not taken.
__m_mform_matrix_MOD___final_m_mform_matrix_Massmatrix:
✗ Branch 3 → 4 not taken.
✗ Branch 3 → 7 not taken.
✗ Branch 4 → 5 not taken.
✗ Branch 4 → 6 not taken.
✗ Branch 8 → 9 not taken.
✗ Branch 8 → 17 not taken.
✗ Branch 10 → 11 not taken.
✗ Branch 10 → 12 not taken.
✗ Branch 12 → 13 not taken.
✗ Branch 12 → 16 not taken.
✗ Branch 13 → 14 not taken.
✗ Branch 13 → 15 not taken.
__m_mform_matrix_MOD___final_m_mform_matrix_Zeromatrix:
✗ Branch 3 → 4 not taken.
✓ Branch 3 → 7 taken 24 times.
✗ Branch 4 → 5 not taken.
✗ Branch 4 → 6 not taken.
✓ Branch 8 → 9 taken 24 times.
✓ Branch 8 → 17 taken 24 times.
✗ Branch 10 → 11 not taken.
✓ Branch 10 → 12 taken 24 times.
✓ Branch 12 → 13 taken 24 times.
✗ Branch 12 → 16 not taken.
✓ Branch 13 → 14 taken 16 times.
✓ Branch 13 → 15 taken 8 times.
|
642 | end module m_mform_matrix |
| 885 |