domain/m_domain_mform.F90
| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | !> @brief Module for determining global indices of a vector of distributed tensor product B-spline spaces | ||
| 2 | module m_domain_mform | ||
| 3 | implicit none | ||
| 4 | |||
| 5 | private | ||
| 6 | public :: determine_linear_indices | ||
| 7 | |||
| 8 | contains | ||
| 9 | |||
| 10 | !> @brief Determine the linear indices for a vector of distributed tensor product B-spline spaces | ||
| 11 | !> | ||
| 12 | !> @param[out] rank_l0 The starting linear index for the combined space | ||
| 13 | !> @param[out] rank_l1 The ending linear index for the combined space | ||
| 14 | !> @param[inout] spaces The vector of TensorProdSpace objects for which the linear indices are determined | ||
| 15 |
1/2✗ Branch 2 → 3 not taken.
✓ Branch 2 → 4 taken 27764 times.
|
27764 | subroutine determine_linear_indices(rank_l0, rank_l1, spaces) |
| 16 | use m_tensorprod_basis, only: TensorProdSpace | ||
| 17 | use m_domain_decomp, only: TensorProdDomain, resp_interval_indices_1d, resp_bspline_indices_1d | ||
| 18 | use m_tensorprod_indices, only: TensorProdIndices, size | ||
| 19 | implicit none | ||
| 20 | |||
| 21 | integer, intent(out) :: rank_l0, rank_l1 | ||
| 22 | type(TensorProdSpace), intent(inout) :: spaces(:) | ||
| 23 | |||
| 24 | integer :: blk, rsi, rsj, rsk | ||
| 25 | |||
| 26 |
2/2✓ Branch 4 → 5 taken 13869 times.
✓ Branch 4 → 6 taken 13895 times.
|
27764 | if (size(spaces) == 1) then |
| 27 | 13869 | rank_l0 = spaces(1)%rank_l0 | |
| 28 | 13869 | rank_l1 = spaces(1)%rank_l1 | |
| 29 | return | ||
| 30 | end if | ||
| 31 | |||
| 32 | ! The unknowns are linearly ordered as follows: | ||
| 33 | ! i, j, k, blk, rank, | ||
| 34 | ! where blk = 1:size(spaces) and we note that the rank is contiguous per shared memory node | ||
| 35 | ! print*,'rank:', spaces(1)%domain%my_rank_i, spaces(1)%domain%my_subinterval_ijk(2), spaces(1)%domain%my_subinterval_ijk(3) | ||
| 36 | |||
| 37 | ! print*,'Initial linear indices:' | ||
| 38 | ! print*,' blk 1: rank_l0:', spaces(1)%rank_l0, ' rank_l1:', spaces(1)%rank_l1 | ||
| 39 | ! print*,' blk 2: rank_l0:', spaces(2)%rank_l0, ' rank_l1:', spaces(2)%rank_l1 | ||
| 40 | ! print*,' blk 3: rank_l0:', spaces(3)%rank_l0, ' rank_l1:', spaces(3)%rank_l1 | ||
| 41 |
2/2✓ Branch 6 → 7 taken 27790 times.
✓ Branch 6 → 20 taken 13895 times.
|
41685 | do blk = 2, size(spaces) |
| 42 | 27790 | spaces(1)%rank_l0 = spaces(1)%rank_l0 + spaces(blk)%rank_l0 | |
| 43 |
2/2✓ Branch 8 → 9 taken 93326 times.
✓ Branch 8 → 18 taken 27790 times.
|
135011 | do rsk = -spaces(1)%domain%nr_eff_neighbours(3), spaces(1)%domain%nr_eff_neighbours(3) |
| 44 |
2/2✓ Branch 10 → 11 taken 189582 times.
✓ Branch 10 → 16 taken 93326 times.
|
310698 | do rsj = -spaces(1)%domain%nr_eff_neighbours(2), spaces(1)%domain%nr_eff_neighbours(2) |
| 45 |
2/2✓ Branch 12 → 13 taken 380046 times.
✓ Branch 12 → 14 taken 189582 times.
|
662954 | do rsi = -spaces(1)%domain%nr_eff_neighbours(1), spaces(1)%domain%nr_eff_neighbours(1) |
| 46 | spaces(1)%neighbours(rsi, rsj, rsk)%their_l0 = & | ||
| 47 | 569628 | spaces(1)%neighbours(rsi, rsj, rsk)%their_l0 + spaces(blk)%neighbours(rsi, rsj, rsk)%their_l0 | |
| 48 | end do | ||
| 49 | end do | ||
| 50 | end do | ||
| 51 | end do | ||
| 52 | 13895 | spaces(1)%rank_l1 = spaces(1)%rank_l0 + size(spaces(1)%rank_resp_bspline) - 1 | |
| 53 |
2/2✓ Branch 22 → 23 taken 46663 times.
✓ Branch 22 → 32 taken 13895 times.
|
60558 | do rsk = -spaces(1)%domain%nr_eff_neighbours(3), spaces(1)%domain%nr_eff_neighbours(3) |
| 54 |
2/2✓ Branch 24 → 25 taken 94791 times.
✓ Branch 24 → 30 taken 46663 times.
|
155349 | do rsj = -spaces(1)%domain%nr_eff_neighbours(2), spaces(1)%domain%nr_eff_neighbours(2) |
| 55 |
2/2✓ Branch 26 → 27 taken 190023 times.
✓ Branch 26 → 28 taken 94791 times.
|
331477 | do rsi = -spaces(1)%domain%nr_eff_neighbours(1), spaces(1)%domain%nr_eff_neighbours(1) |
| 56 | spaces(1)%neighbours(rsi, rsj, rsk)%their_l1 = & | ||
| 57 | 284814 | spaces(1)%neighbours(rsi, rsj, rsk)%their_l0 + size(spaces(1)%neighbours(rsi, rsj, rsk)%their_rank_resp_bspline) - 1 | |
| 58 | end do | ||
| 59 | end do | ||
| 60 | end do | ||
| 61 | |||
| 62 |
2/2✓ Branch 33 → 34 taken 27790 times.
✓ Branch 33 → 47 taken 13895 times.
|
41685 | do blk = 2, size(spaces) |
| 63 | 27790 | spaces(blk)%rank_l0 = spaces(blk - 1)%rank_l1 + 1 | |
| 64 | 27790 | spaces(blk)%rank_l1 = spaces(blk)%rank_l0 + size(spaces(blk)%rank_resp_bspline) - 1 | |
| 65 | |||
| 66 |
2/2✓ Branch 35 → 36 taken 93326 times.
✓ Branch 35 → 45 taken 27790 times.
|
135011 | do rsk = -spaces(1)%domain%nr_eff_neighbours(3), spaces(1)%domain%nr_eff_neighbours(3) |
| 67 |
2/2✓ Branch 37 → 38 taken 189582 times.
✓ Branch 37 → 43 taken 93326 times.
|
310698 | do rsj = -spaces(1)%domain%nr_eff_neighbours(2), spaces(1)%domain%nr_eff_neighbours(2) |
| 68 |
2/2✓ Branch 39 → 40 taken 380046 times.
✓ Branch 39 → 41 taken 189582 times.
|
662954 | do rsi = -spaces(1)%domain%nr_eff_neighbours(1), spaces(1)%domain%nr_eff_neighbours(1) |
| 69 | spaces(blk)%neighbours(rsi, rsj, rsk)%their_l0 = & | ||
| 70 | 380046 | spaces(blk - 1)%neighbours(rsi, rsj, rsk)%their_l1 + 1 | |
| 71 | spaces(blk)%neighbours(rsi, rsj, rsk)%their_l1 = & | ||
| 72 | 569628 | spaces(blk)%neighbours(rsi, rsj, rsk)%their_l0 + size(spaces(blk)%neighbours(rsi, rsj, rsk)%their_rank_resp_bspline) - 1 | |
| 73 | end do | ||
| 74 | end do | ||
| 75 | end do | ||
| 76 | end do | ||
| 77 | |||
| 78 | 13895 | rank_l0 = spaces(1)%rank_l0 | |
| 79 | 13895 | rank_l1 = spaces(size(spaces))%rank_l1 | |
| 80 | |||
| 81 | ! Set empty blocks to 0, -1 | ||
| 82 |
2/2✓ Branch 49 → 50 taken 41685 times.
✓ Branch 49 → 67 taken 13895 times.
|
55580 | do blk = 1, size(spaces) |
| 83 |
2/2✓ Branch 50 → 51 taken 4719 times.
✓ Branch 50 → 52 taken 36966 times.
|
41685 | if (spaces(blk)%rank_l0 > spaces(blk)%rank_l1) then |
| 84 | 4719 | spaces(blk)%rank_l0 = 0 | |
| 85 | 4719 | spaces(blk)%rank_l1 = -1 | |
| 86 | end if | ||
| 87 | |||
| 88 |
2/2✓ Branch 53 → 54 taken 139989 times.
✓ Branch 53 → 65 taken 41685 times.
|
195569 | do rsk = -spaces(1)%domain%nr_eff_neighbours(3), spaces(1)%domain%nr_eff_neighbours(3) |
| 89 |
2/2✓ Branch 55 → 56 taken 284373 times.
✓ Branch 55 → 63 taken 139989 times.
|
466047 | do rsj = -spaces(1)%domain%nr_eff_neighbours(2), spaces(1)%domain%nr_eff_neighbours(2) |
| 90 |
2/2✓ Branch 57 → 58 taken 570069 times.
✓ Branch 57 → 61 taken 284373 times.
|
994431 | do rsi = -spaces(1)%domain%nr_eff_neighbours(1), spaces(1)%domain%nr_eff_neighbours(1) |
| 91 |
2/2✓ Branch 58 → 59 taken 566639 times.
✓ Branch 58 → 60 taken 3430 times.
|
854442 | if (spaces(blk)%neighbours(rsi, rsj, rsk)%their_l0 > spaces(blk)%neighbours(rsi, rsj, rsk)%their_l1) then |
| 92 | 566639 | spaces(blk)%neighbours(rsi, rsj, rsk)%their_l0 = 0 | |
| 93 | 566639 | spaces(blk)%neighbours(rsi, rsj, rsk)%their_l1 = -1 | |
| 94 | end if | ||
| 95 | end do | ||
| 96 | end do | ||
| 97 | end do | ||
| 98 | end do | ||
| 99 | |||
| 100 | ! print*,'Updated linear indices:' | ||
| 101 | ! print*,' blk 1: rank_l0:', spaces(1)%rank_l0, ' rank_l1:', spaces(1)%rank_l1 | ||
| 102 | ! print*,' blk 2: rank_l0:', spaces(2)%rank_l0, ' rank_l1:', spaces(2)%rank_l1 | ||
| 103 | ! print*,' blk 3: rank_l0:', spaces(3)%rank_l0, ' rank_l1:', spaces(3)%rank_l1 | ||
| 104 | |||
| 105 | end subroutine determine_linear_indices | ||
| 106 | |||
| 107 | end module m_domain_mform | ||
| 108 |