FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_solver_las_11.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 
7 
8  private
9 
10  public :: hecmw_matvec_11
11  public :: hecmw_matresid_11
12  public :: hecmw_rel_resid_l2_11
13 
14 contains
15 
16  !C
17  !C***
18  !C*** hecmw_matvec_11
19  !C***
20  !C
21  subroutine hecmw_matvec_11 (hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
23  use m_hecmw_comm_f
24 
25  implicit none
26  type (hecmwst_local_mesh), intent(in) :: hecmesh
27  type (hecmwst_matrix), intent(in) :: hecmat
28  real(kind=kreal), intent(in) :: x(:)
29  real(kind=kreal), intent(out) :: y(:)
30  real(kind=kreal), intent(inout) :: time_ax
31  real(kind=kreal), intent(inout), optional :: commtime
32 
33  real(kind=kreal) :: start_time, end_time
34  integer(kind=kint) :: i, j, js, je, in
35  real(kind=kreal) :: yv
36 
37  start_time= hecmw_wtime()
38  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
39  end_time= hecmw_wtime()
40  if (present(commtime)) commtime = commtime + end_time - start_time
41 
42  start_time= hecmw_wtime()
43  do i= 1, hecmat%N
44  yv= hecmat%D(i) * x(i)
45  js= hecmat%indexL(i-1) + 1
46  je= hecmat%indexL(i )
47  do j= js, je
48  in= hecmat%itemL(j)
49  yv= yv + hecmat%AL(j) * x(in)
50  enddo
51  js= hecmat%indexU(i-1) + 1
52  je= hecmat%indexU(i )
53  do j= js, je
54  in= hecmat%itemU(j)
55  yv= yv + hecmat%AU(j) * x(in)
56  enddo
57  y(i)= yv
58  enddo
59  end_time = hecmw_wtime()
60  time_ax = time_ax + end_time - start_time
61 
62  end subroutine hecmw_matvec_11
63 
64  !C
65  !C***
66  !C*** hecmw_matresid_11
67  !C***
68  !C
69  subroutine hecmw_matresid_11 (hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
71 
72  implicit none
73  real(kind=kreal) :: x(:), b(:), r(:)
74  type (hecmwst_matrix) :: hecmat
75  type (hecmwst_local_mesh) :: hecmesh
76  real(kind=kreal) :: time_ax
77  real(kind=kreal), optional :: commtime
78 
79  integer(kind=kint) :: i
80  real(kind=kreal) :: tcomm
81 
82  tcomm = 0.d0
83  call hecmw_matvec_11 (hecmesh, hecmat, x, r, time_ax, tcomm)
84  if (present(commtime)) commtime = commtime + tcomm
85 
86  do i = 1, hecmat%N
87  r(i) = b(i) - r(i)
88  enddo
89 
90  end subroutine hecmw_matresid_11
91 
92  !C
93  !C***
94  !C*** hecmw_rel_resid_L2_11
95  !C***
96  !C
97  function hecmw_rel_resid_l2_11 (hecMESH, hecMAT, time_Ax, COMMtime)
100 
101  implicit none
102  real(kind=kreal) :: hecmw_rel_resid_l2_11
103  type ( hecmwst_local_mesh ), intent(in) :: hecmesh
104  type ( hecmwst_matrix ), intent(in) :: hecmat
105  real(kind=kreal) :: time_ax
106  real(kind=kreal), optional :: commtime
107 
108  real(kind=kreal) :: r(hecmat%NDOF*hecmat%NP)
109  real(kind=kreal) :: bnorm2, rnorm2
110  real(kind=kreal) :: tcomm
111 
112  tcomm = 0.d0
113  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, hecmat%B, hecmat%B, bnorm2, tcomm)
114  if (bnorm2 == 0.d0) then
115  bnorm2 = 1.d0
116  endif
117  call hecmw_matresid_11(hecmesh, hecmat, hecmat%X, hecmat%B, r, time_ax, tcomm)
118  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, r, r, rnorm2, tcomm)
119  if (present(commtime)) commtime = commtime + tcomm
120  hecmw_rel_resid_l2_11 = sqrt(rnorm2 / bnorm2)
121 
122  end function hecmw_rel_resid_l2_11
123 
124 end module hecmw_solver_las_11
hecmw_solver_las_11::hecmw_matvec_11
subroutine, public hecmw_matvec_11(hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
Definition: hecmw_solver_las_11.f90:22
hecmw_util::hecmw_wtime
real(kind=kreal) function hecmw_wtime()
Definition: hecmw_util_f.F90:549
hecmw_solver_misc::hecmw_innerproduct_r
subroutine hecmw_innerproduct_r(hecMESH, ndof, X, Y, sum, COMMtime)
Definition: hecmw_solver_misc.f90:47
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_util::hecmwst_local_mesh
Definition: hecmw_util_f.F90:234
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_solver_las_11::hecmw_rel_resid_l2_11
real(kind=kreal) function, public hecmw_rel_resid_l2_11(hecMESH, hecMAT, time_Ax, COMMtime)
Definition: hecmw_solver_las_11.f90:98
hecmw_solver_las_11
Definition: hecmw_solver_las_11.f90:6
hecmw_solver_misc
Definition: hecmw_solver_misc.f90:6
hecmw_solver_las_11::hecmw_matresid_11
subroutine, public hecmw_matresid_11(hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
Definition: hecmw_solver_las_11.f90:70
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444