|
| 1 | +!>@file sph_zero_degree_matrices.f90 |
| 2 | +!!@brief module sph_zero_degree_matrices |
| 3 | +!! |
| 4 | +!!@author H. Matsui |
| 5 | +!!@date Programmed in Apr, 2009 |
| 6 | +! |
| 7 | +!>@brief Construct matrix for degree 0 |
| 8 | +!! |
| 9 | +!!@verbatim |
| 10 | +!! subroutine copy_to_band3_mat_w_center(nri, c_evo, mat3, mat00_3) |
| 11 | +!! subroutine copy_to_band3_mat_no_center(nri, mat3, mat00_3) |
| 12 | +!! subroutine add_scalar_poisson_mat_fill_ctr(nri, r_CTR1, & |
| 13 | +!! & fdm2_fix_dr_center, fdm2_fix_fld_ctr1, coef_p, mat00_3) |
| 14 | +!! subroutine add_scalar_poisson_mat_fix_ctr(nri, r_CTR1, & |
| 15 | +!! & fdm2_fix_fld_ctr1, coef_p, mat00_3) |
| 16 | +!! subroutine add_scalar_poisson_mat_no_fld(nri, mat00_3) |
| 17 | +!! integer(kind = kint), intent(in) :: nri |
| 18 | +!! real(kind = kreal), intent(in) :: c_evo |
| 19 | +!! real(kind = kreal), intent(in) :: mat3(3,nri) |
| 20 | +!! real(kind = kreal), intent(in) :: r_CTR1(0:2) |
| 21 | +!! real(kind = kreal), intent(in) :: coef_p |
| 22 | +!! real(kind = kreal), intent(in) :: fdm2_fix_dr_center(-1:1,3) |
| 23 | +!! real(kind = kreal), intent(in) :: fdm2_fix_fld_ctr1(-1:1,3) |
| 24 | +!! real(kind = kreal), intent(inout) :: mat00_3(3,0:nri) |
| 25 | +!!@endverbatim |
| 26 | +! |
| 27 | +!!@n @param jmax Number of local spherical harmonics mode |
| 28 | +!!@n @param fdm2_fix_fld_ctr1(-1:1,3) |
| 29 | +!! Matrix to evaluate radial derivative |
| 30 | +!! for center with fixed field |
| 31 | +! |
| 32 | +! |
| 33 | + module sph_zero_degree_matrices |
| 34 | +! |
| 35 | + use m_precision |
| 36 | + use m_constants |
| 37 | +! |
| 38 | + implicit none |
| 39 | +! |
| 40 | +! ----------------------------------------------------------------------- |
| 41 | +! |
| 42 | + contains |
| 43 | +! |
| 44 | +! ----------------------------------------------------------------------- |
| 45 | +! |
| 46 | + subroutine copy_to_band3_mat_w_center(nri, c_evo, mat3, mat00_3) |
| 47 | +! |
| 48 | + integer(kind = kint), intent(in) :: nri |
| 49 | + real(kind = kreal), intent(in) :: c_evo |
| 50 | + real(kind = kreal), intent(in) :: mat3(3,nri) |
| 51 | + real(kind = kreal), intent(inout) :: mat00_3(3,0:nri) |
| 52 | +! |
| 53 | + integer(kind = kint) :: k |
| 54 | +! |
| 55 | +! |
| 56 | + if(c_evo .eq. zero) then |
| 57 | + mat00_3(2,0) = zero |
| 58 | + else |
| 59 | + mat00_3(2,0) = one |
| 60 | + end if |
| 61 | + mat00_3(1,1) = zero |
| 62 | +! |
| 63 | + mat00_3(3,0) = zero |
| 64 | + if(c_evo .eq. zero) then |
| 65 | + mat00_3(2,1) = zero |
| 66 | + else |
| 67 | + mat00_3(2,1) = one |
| 68 | + end if |
| 69 | + mat00_3(1,2) = zero |
| 70 | +! |
| 71 | + do k = 2, nri-1 |
| 72 | + mat00_3(3,k-1) = mat3(3,k-1) |
| 73 | + mat00_3(2,k ) = mat3(2,k ) |
| 74 | + mat00_3(1,k+1) = mat3(1,k+1) |
| 75 | + end do |
| 76 | + mat00_3(3,nri-1) = mat3(3,nri-1) |
| 77 | + mat00_3(2,nri ) = mat3(2,nri ) |
| 78 | +! |
| 79 | + end subroutine copy_to_band3_mat_w_center |
| 80 | +! |
| 81 | +! ----------------------------------------------------------------------- |
| 82 | +! |
| 83 | + subroutine copy_to_band3_mat_no_center(nri, mat3, mat00_3) |
| 84 | +! |
| 85 | + integer(kind = kint), intent(in) :: nri |
| 86 | + real(kind = kreal), intent(in) :: mat3(3,nri) |
| 87 | + real(kind = kreal), intent(inout) :: mat00_3(3,0:nri) |
| 88 | +! |
| 89 | + integer(kind = kint) :: k |
| 90 | +! |
| 91 | +! |
| 92 | + mat00_3(2,0) = one |
| 93 | + mat00_3(1,1) = zero |
| 94 | +! |
| 95 | + mat00_3(3,0) = zero |
| 96 | + mat00_3(2,1) = mat3(2,1) |
| 97 | + mat00_3(1,2) = mat3(1,2) |
| 98 | +! |
| 99 | + do k = 2, nri-1 |
| 100 | + mat00_3(3,k-1) = mat3(3,k-1) |
| 101 | + mat00_3(2,k ) = mat3(2,k ) |
| 102 | + mat00_3(1,k+1) = mat3(1,k+1) |
| 103 | + end do |
| 104 | + mat00_3(3,nri-1) = mat3(3,nri-1) |
| 105 | + mat00_3(2,nri ) = mat3(2,nri ) |
| 106 | +! |
| 107 | + end subroutine copy_to_band3_mat_no_center |
| 108 | +! |
| 109 | +! ----------------------------------------------------------------------- |
| 110 | +! |
| 111 | + subroutine add_scalar_poisson_mat_fill_ctr(nri, r_CTR1, & |
| 112 | + & fdm2_fix_dr_center, fdm2_fix_fld_ctr1, coef_p, mat00_3) |
| 113 | +! |
| 114 | + integer(kind = kint), intent(in) :: nri |
| 115 | + real(kind = kreal), intent(in) :: r_CTR1(0:2) |
| 116 | + real(kind = kreal), intent(in) :: coef_p |
| 117 | + real(kind = kreal), intent(in) :: fdm2_fix_dr_center(-1:1,3) |
| 118 | + real(kind = kreal), intent(in) :: fdm2_fix_fld_ctr1(-1:1,3) |
| 119 | +! |
| 120 | + real(kind = kreal), intent(inout) :: mat00_3(3,0:nri) |
| 121 | +! |
| 122 | +! |
| 123 | +! |
| 124 | +! mat00_3(3,-1) = mat00_3(3,-1) - coef_p * fdm2_fix_dr_center(-1,3) |
| 125 | + mat00_3(2,0) = mat00_3(2,0) - coef_p * fdm2_fix_dr_center(0,3) |
| 126 | + mat00_3(1,1) = mat00_3(1,1) - coef_p * fdm2_fix_dr_center(1,3) |
| 127 | +! |
| 128 | + mat00_3(3,0) = mat00_3(3,0) - coef_p * (fdm2_fix_fld_ctr1(-1,3) & |
| 129 | + & + two*r_CTR1(1) * fdm2_fix_fld_ctr1(-1,2)) |
| 130 | + mat00_3(2,1) = mat00_3(2,1) - coef_p * (fdm2_fix_fld_ctr1(0,3) & |
| 131 | + & + two*r_CTR1(1) * fdm2_fix_fld_ctr1(0,2)) |
| 132 | + mat00_3(1,2) = mat00_3(1,2) - coef_p * (fdm2_fix_fld_ctr1(1,3) & |
| 133 | + & + two*r_CTR1(1) * fdm2_fix_fld_ctr1(1,2)) |
| 134 | +! |
| 135 | + end subroutine add_scalar_poisson_mat_fill_ctr |
| 136 | +! |
| 137 | +! ----------------------------------------------------------------------- |
| 138 | +! |
| 139 | + subroutine add_scalar_poisson_mat_fix_ctr(nri, r_CTR1, & |
| 140 | + & fdm2_fix_fld_ctr1, coef_p, mat00_3) |
| 141 | +! |
| 142 | + integer(kind = kint), intent(in) :: nri |
| 143 | + real(kind = kreal), intent(in) :: r_CTR1(0:2) |
| 144 | + real(kind = kreal), intent(in) :: coef_p |
| 145 | + real(kind = kreal), intent(in) :: fdm2_fix_fld_ctr1(-1:1,3) |
| 146 | +! |
| 147 | + real(kind = kreal), intent(inout) :: mat00_3(3,0:nri) |
| 148 | +! |
| 149 | +! |
| 150 | + mat00_3(2,0) = one |
| 151 | + mat00_3(1,1) = zero |
| 152 | +! |
| 153 | + mat00_3(3,0) = mat00_3(3,0) - coef_p * (fdm2_fix_fld_ctr1(-1,3) & |
| 154 | + & + two*r_CTR1(1) * fdm2_fix_fld_ctr1(-1,2)) |
| 155 | + mat00_3(2,1) = mat00_3(2,1) - coef_p * (fdm2_fix_fld_ctr1(0,3) & |
| 156 | + & + two*r_CTR1(1) * fdm2_fix_fld_ctr1(0,2)) |
| 157 | + mat00_3(1,2) = mat00_3(1,2) - coef_p * (fdm2_fix_fld_ctr1(1,3) & |
| 158 | + & + two*r_CTR1(1) * fdm2_fix_fld_ctr1(1,2)) |
| 159 | +! |
| 160 | + end subroutine add_scalar_poisson_mat_fix_ctr |
| 161 | +! |
| 162 | +! ----------------------------------------------------------------------- |
| 163 | +! |
| 164 | + subroutine add_scalar_poisson_mat_no_fld(nri, mat00_3) |
| 165 | +! |
| 166 | + integer(kind = kint), intent(in) :: nri |
| 167 | + real(kind = kreal), intent(inout) :: mat00_3(3,0:nri) |
| 168 | +! |
| 169 | +! |
| 170 | + mat00_3(2,0) = one |
| 171 | + mat00_3(1,1) = zero |
| 172 | +! |
| 173 | + mat00_3(3,0) = zero |
| 174 | + mat00_3(2,1) = one |
| 175 | + mat00_3(1,2) = zero |
| 176 | +! |
| 177 | + end subroutine add_scalar_poisson_mat_no_fld |
| 178 | +! |
| 179 | +! ----------------------------------------------------------------------- |
| 180 | +! |
| 181 | + end module sph_zero_degree_matrices |
0 commit comments