FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_ML_helper_nn_f.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
5 
6 subroutine hecmw_ml_getrow_nn(id, n_requested_rows, requested_rows, &
7  allocated_space, cols, values, row_lengths, ierr)
9  use hecmw_mat_id
10  implicit none
11  integer(kind=kint), intent(in) :: id
12  integer(kind=kint), intent(in) :: n_requested_rows
13  integer(kind=kint), intent(in) :: requested_rows(n_requested_rows)
14  integer(kind=kint), intent(in) :: allocated_space
15  integer(kind=kint), intent(out) :: cols(allocated_space)
16  real(kind=kreal), intent(out) :: values(allocated_space)
17  integer(kind=kint), intent(out) :: row_lengths(n_requested_rows)
18  integer(kind=kint), intent(out) :: ierr
19  type(hecmwst_matrix), pointer :: hecMAT
20  type(hecmwst_local_mesh), pointer :: hecMESH
21  integer(kind=kint) :: m, i, row, inod, idof, nl, nd, nu, js, je, j, jj, jdof, start, ndof
22  ierr = 0
23  call hecmw_mat_id_get(id, hecmat, hecmesh)
24  ndof = hecmat%NDOF
25  m = 1
26  do i = 1, n_requested_rows
27  row = requested_rows(i) + 1 ! '+1' for Fortran-numbering
28  inod = (row-1)/ndof + 1
29  idof = row - (inod-1)*ndof
30  nl = (hecmat%indexL(inod) - hecmat%indexL(inod-1)) * ndof
31  nd = ndof
32  nu = (hecmat%indexU(inod) - hecmat%indexU(inod-1)) * ndof
33  if (allocated_space < m + nl + nd + nu) return
34  start = m
35  js = hecmat%indexL(inod-1)+1
36  je = hecmat%indexL(inod)
37  do j = js, je
38  jj = hecmat%itemL(j)
39  do jdof = 1, ndof
40  cols(m) = (jj-1)*ndof + jdof - 1 ! '-1' for C-numbering
41  values(m) = hecmat%AL((j-1)*ndof*ndof + (idof-1)*ndof + jdof)
42  m = m+1
43  enddo
44  enddo
45  do jdof = 1, ndof
46  cols(m) = (inod-1)*ndof + jdof - 1 ! '-1' for C-numbering
47  values(m) = hecmat%D((inod-1)*ndof*ndof + (idof-1)*ndof + jdof)
48  m = m+1
49  enddo
50  js = hecmat%indexU(inod-1)+1
51  je = hecmat%indexU(inod)
52  do j = js, je
53  jj = hecmat%itemU(j)
54  do jdof = 1, ndof
55  cols(m) = (jj-1)*ndof + jdof - 1 ! '-1' for C-numbering
56  values(m) = hecmat%AU((j-1)*ndof*ndof + (idof-1)*ndof + jdof)
57  m = m+1
58  enddo
59  enddo
60  row_lengths(i) = m - start
61  enddo
62  ierr = 1
63 end subroutine hecmw_ml_getrow_nn
64 
65 subroutine hecmw_ml_matvec_nn(id, in_length, p, out_length, ap, ierr)
67  use hecmw_mat_id
69  implicit none
70  integer(kind=kint), intent(in) :: id
71  integer(kind=kint), intent(in) :: in_length
72  real(kind=kreal), intent(in) :: p(in_length)
73  integer(kind=kint), intent(in) :: out_length
74  real(kind=kreal), intent(out) :: ap(out_length)
75  integer(kind=kint), intent(out) :: ierr
76  type(hecmwst_matrix), pointer :: hecMAT
77  type(hecmwst_local_mesh), pointer :: hecMESH
78  real(kind=kreal), allocatable :: w(:)
79  integer(kind=kint) :: i
80  call hecmw_mat_id_get(id, hecmat, hecmesh)
81  allocate(w(hecmat%NP*hecmat%NDOF))
82  do i = 1, hecmat%N*hecmat%NDOF
83  w(i) = p(i)
84  enddo
85  call hecmw_matvec(hecmesh, hecmat, w, ap)
86  deallocate(w)
87  ierr = 0
88 end subroutine hecmw_ml_matvec_nn
89 
90 subroutine hecmw_ml_comm_nn(id, x, ierr)
92  use hecmw_mat_id
93  use m_hecmw_comm_f
94  implicit none
95  integer(kind=kint), intent(in) :: id
96  real(kind=kreal), intent(inout) :: x(*)
97  integer(kind=kint), intent(out) :: ierr
98  type(hecmwst_matrix), pointer :: hecMAT
99  type(hecmwst_local_mesh), pointer :: hecMESH
100  call hecmw_mat_id_get(id, hecmat, hecmesh)
101  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
102  ierr = 0
103 end subroutine hecmw_ml_comm_nn
104 
105 subroutine hecmw_ml_smoother_diag_setup_nn(id, ierr)
107  use hecmw_mat_id
109  implicit none
110  integer(kind=kint), intent(in) :: id
111  integer(kind=kint), intent(out) :: ierr
112  type(hecmwst_matrix), pointer :: hecMAT
113  type(hecmwst_local_mesh), pointer :: hecMESH
114  call hecmw_mat_id_get(id, hecmat, hecmesh)
115  call hecmw_precond_diag_nn_setup(hecmat)
116  ierr = 0
117 end subroutine hecmw_ml_smoother_diag_setup_nn
118 
119 subroutine hecmw_ml_smoother_diag_apply_nn(id, x_length, x, rhs_length, rhs, ierr)
121  use hecmw_mat_id
123  use hecmw_solver_las
125  implicit none
126  integer(kind=kint), intent(in) :: id
127  integer(kind=kint), intent(in) :: x_length
128  real(kind=kreal), intent(inout) :: x(x_length)
129  integer(kind=kint), intent(in) :: rhs_length
130  real(kind=kreal), intent(in) :: rhs(rhs_length)
131  integer(kind=kint), intent(out) :: ierr
132  type(hecmwst_matrix), pointer :: hecMAT
133  type(hecmwst_local_mesh), pointer :: hecMESH
134 
135  real(kind=kreal), allocatable :: resid(:)
136  integer(kind=kint) :: i
137  real(kind=kreal) :: commtime
138  integer(kind=kint) :: num_sweeps, i_sweep
139  integer(kind=kint) :: opt(10)
140 
141  call hecmw_mat_id_get(id, hecmat, hecmesh)
142  call hecmw_mat_get_solver_opt(hecmat, opt)
143  num_sweeps = opt(6)
144  allocate(resid(hecmat%NP * hecmat%NDOF))
145  do i_sweep = 1, num_sweeps
146  ! {resid} = {rhs} - [A] {x}
147  call hecmw_matresid(hecmesh, hecmat, x, rhs, resid, commtime)
148  ! {delta_x} = [M]^-1 {resid}
149  call hecmw_precond_diag_nn_apply(resid, hecmat%NDOF)
150  ! {x} = {x} + {delta_x}
151  do i=1,x_length
152  x(i) = x(i) + resid(i)
153  enddo
154  enddo
155  deallocate(resid)
156  ierr = 0
157 end subroutine hecmw_ml_smoother_diag_apply_nn
158 
159 subroutine hecmw_ml_smoother_diag_clear_nn(id, ierr)
161  use hecmw_mat_id
163  implicit none
164  integer(kind=kint), intent(in) :: id
165  integer(kind=kint), intent(out) :: ierr
166  type(hecmwst_matrix), pointer :: hecMAT
167  type(hecmwst_local_mesh), pointer :: hecMESH
168  call hecmw_mat_id_get(id, hecmat, hecmesh)
170  ierr = 0
171 end subroutine hecmw_ml_smoother_diag_clear_nn
172 
173 subroutine hecmw_ml_smoother_ssor_setup_nn(id, ierr)
175  use hecmw_mat_id
177  implicit none
178  integer(kind=kint), intent(in) :: id
179  integer(kind=kint), intent(out) :: ierr
180  type(hecmwst_matrix), pointer :: hecMAT
181  type(hecmwst_local_mesh), pointer :: hecMESH
182  call hecmw_mat_id_get(id, hecmat, hecmesh)
183  call hecmw_precond_ssor_nn_setup(hecmat)
184  ierr = 0
185 end subroutine hecmw_ml_smoother_ssor_setup_nn
186 
187 subroutine hecmw_ml_smoother_ssor_apply_nn(id, x_length, x, rhs_length, rhs, ierr)
189  use hecmw_mat_id
191  use hecmw_solver_las
193  implicit none
194  integer(kind=kint), intent(in) :: id
195  integer(kind=kint), intent(in) :: x_length
196  real(kind=kreal), intent(inout) :: x(x_length)
197  integer(kind=kint), intent(in) :: rhs_length
198  real(kind=kreal), intent(in) :: rhs(rhs_length)
199  integer(kind=kint), intent(out) :: ierr
200  type(hecmwst_matrix), pointer :: hecMAT
201  type(hecmwst_local_mesh), pointer :: hecMESH
202 
203  real(kind=kreal), allocatable :: resid(:)
204  integer(kind=kint) :: i
205  real(kind=kreal) :: commtime
206  integer(kind=kint) :: num_sweeps, i_sweep
207  integer(kind=kint) :: opt(10)
208 
209  call hecmw_mat_id_get(id, hecmat, hecmesh)
210  call hecmw_mat_get_solver_opt(hecmat, opt)
211  num_sweeps = opt(6)
212  allocate(resid(hecmat%NP * hecmat%NDOF))
213  do i_sweep = 1, num_sweeps
214  ! {resid} = {rhs} - [A] {x}
215  call hecmw_matresid(hecmesh, hecmat, x, rhs, resid, commtime)
216  ! {delta_x} = [M]^-1 {resid}
217  call hecmw_precond_ssor_nn_apply(resid, hecmat%NDOF)
218  ! {x} = {x} + {delta_x}
219  do i=1,x_length
220  x(i) = x(i) + resid(i)
221  enddo
222  enddo
223  deallocate(resid)
224  ierr = 0
225 end subroutine hecmw_ml_smoother_ssor_apply_nn
226 
227 subroutine hecmw_ml_smoother_ssor_clear_nn(id, ierr)
229  use hecmw_mat_id
231  implicit none
232  integer(kind=kint), intent(in) :: id
233  integer(kind=kint), intent(out) :: ierr
234  type(hecmwst_matrix), pointer :: hecMAT
235  type(hecmwst_local_mesh), pointer :: hecMESH
236  call hecmw_mat_id_get(id, hecmat, hecmesh)
237  call hecmw_precond_ssor_nn_clear(hecmat)
238  ierr = 0
239 end subroutine hecmw_ml_smoother_ssor_clear_nn
hecmw_precond_ssor_nn::hecmw_precond_ssor_nn_setup
subroutine, public hecmw_precond_ssor_nn_setup(hecMAT)
Definition: hecmw_precond_SSOR_nn.f90:47
hecmw_precond_diag_nn::hecmw_precond_diag_nn_clear
subroutine, public hecmw_precond_diag_nn_clear()
Definition: hecmw_precond_DIAG_nn.f90:131
hecmw_ml_matvec_nn
subroutine hecmw_ml_matvec_nn(id, in_length, p, out_length, ap, ierr)
Definition: hecmw_ML_helper_nn_f.f90:66
hecmw_solver_las::hecmw_matvec
subroutine, public hecmw_matvec(hecMESH, hecMAT, X, Y, COMMtime)
Definition: hecmw_solver_las.f90:43
hecmw_matrix_misc::hecmw_mat_get_solver_opt
subroutine, public hecmw_mat_get_solver_opt(hecMAT, solver_opt)
Definition: hecmw_matrix_misc.f90:657
hecmw_precond_diag_nn
Definition: hecmw_precond_DIAG_nn.f90:11
hecmw_precond_diag_nn::hecmw_precond_diag_nn_setup
subroutine, public hecmw_precond_diag_nn_setup(hecMAT)
Definition: hecmw_precond_DIAG_nn.f90:28
hecmw_solver_las
Definition: hecmw_solver_las.f90:6
hecmw_ml_smoother_ssor_setup_nn
subroutine hecmw_ml_smoother_ssor_setup_nn(id, ierr)
Definition: hecmw_ML_helper_nn_f.f90:174
hecmw_solver_las::hecmw_matresid
subroutine, public hecmw_matresid(hecMESH, hecMAT, X, B, R, COMMtime)
Definition: hecmw_solver_las.f90:95
hecmw_ml_smoother_ssor_apply_nn
subroutine hecmw_ml_smoother_ssor_apply_nn(id, x_length, x, rhs_length, rhs, ierr)
Definition: hecmw_ML_helper_nn_f.f90:188
hecmw_ml_getrow_nn
subroutine hecmw_ml_getrow_nn(id, n_requested_rows, requested_rows, allocated_space, cols, values, row_lengths, ierr)
Definition: hecmw_ML_helper_nn_f.f90:8
hecmw_precond_diag_nn::hecmw_precond_diag_nn_apply
subroutine, public hecmw_precond_diag_nn_apply(WW, NDOF)
Definition: hecmw_precond_DIAG_nn.f90:96
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_util::hecmwst_local_mesh
Definition: hecmw_util_f.F90:234
hecmw_ml_smoother_ssor_clear_nn
subroutine hecmw_ml_smoother_ssor_clear_nn(id, ierr)
Definition: hecmw_ML_helper_nn_f.f90:228
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_precond_ssor_nn
Definition: hecmw_precond_SSOR_nn.f90:11
hecmw_ml_smoother_diag_apply_nn
subroutine hecmw_ml_smoother_diag_apply_nn(id, x_length, x, rhs_length, rhs, ierr)
Definition: hecmw_ML_helper_nn_f.f90:120
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_mat_id
Definition: hecmw_mat_id.f90:6
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_ml_smoother_diag_clear_nn
subroutine hecmw_ml_smoother_diag_clear_nn(id, ierr)
Definition: hecmw_ML_helper_nn_f.f90:160
hecmw_ml_smoother_diag_setup_nn
subroutine hecmw_ml_smoother_diag_setup_nn(id, ierr)
Definition: hecmw_ML_helper_nn_f.f90:106
hecmw_precond_ssor_nn::hecmw_precond_ssor_nn_apply
subroutine, public hecmw_precond_ssor_nn_apply(ZP, NDOF)
Definition: hecmw_precond_SSOR_nn.f90:172
hecmw_precond_ssor_nn::hecmw_precond_ssor_nn_clear
subroutine, public hecmw_precond_ssor_nn_clear(hecMAT)
Definition: hecmw_precond_SSOR_nn.f90:343
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
hecmw_ml_comm_nn
subroutine hecmw_ml_comm_nn(id, x, ierr)
Definition: hecmw_ML_helper_nn_f.f90:91
hecmw_mat_id::hecmw_mat_id_get
subroutine, public hecmw_mat_id_get(id, hecMAT, hecMESH)
Definition: hecmw_mat_id.f90:49
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444