FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_precond_DIAG_22.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 !C
7 !C***
8 !C*** module hecmw_precond_DIAG_22
9 !C***
10 !C
12  use hecmw_util
13 
14  private
15 
19 
20  integer(kind=kint) :: N
21  real(kind=kreal), pointer :: alu(:) => null()
22 
23  logical, save :: INITIALIZED = .false.
24 
25 contains
26 
27  subroutine hecmw_precond_diag_22_setup(hecMAT)
29  implicit none
30  type(hecmwst_matrix), intent(inout) :: hecmat
31  integer(kind=kint ) :: np
32  real (kind=kreal) :: sigma_diag
33  real(kind=kreal), pointer:: d(:)
34 
35  real (kind=kreal):: alutmp(2,2), pw(2)
36  integer(kind=kint ):: ii, i, j, k
37 
38  if (initialized) then
39  if (hecmat%Iarray(98) == 0 .and. hecmat%Iarray(97) == 0) return
41  endif
42 
43  n = hecmat%N
44  np = hecmat%NP
45  d => hecmat%D
46  sigma_diag = hecmw_mat_get_sigma_diag(hecmat)
47 
48  allocate(alu(4*np))
49  alu = 0.d0
50 
51  do ii= 1, n
52  alu(4*ii-3) = d(4*ii-3)
53  alu(4*ii-2) = d(4*ii-2)
54  alu(4*ii-1) = d(4*ii-1)
55  alu(4*ii-0) = d(4*ii-0)
56  enddo
57 
58  !$omp parallel default(none),private(ii,ALUtmp,k,i,j,PW),shared(N,ALU,SIGMA_DIAG)
59  !$omp do
60  do ii= 1, n
61  alutmp(1,1)= alu(4*ii-3) * sigma_diag
62  alutmp(1,2)= alu(4*ii-2)
63  alutmp(2,1)= alu(4*ii-1)
64  alutmp(2,2)= alu(4*ii-0) * sigma_diag
65 
66  do k= 1, 2
67  alutmp(k,k)= 1.d0/alutmp(k,k)
68  do i= k+1, 2
69  alutmp(i,k)= alutmp(i,k) * alutmp(k,k)
70  do j= k+1, 2
71  pw(j)= alutmp(i,j) - alutmp(i,k)*alutmp(k,j)
72  enddo
73  do j= k+1, 2
74  alutmp(i,j)= pw(j)
75  enddo
76  enddo
77  enddo
78  alu(4*ii-3)= alutmp(1,1)
79  alu(4*ii-2)= alutmp(1,2)
80  alu(4*ii-1)= alutmp(2,1)
81  alu(4*ii-0)= alutmp(2,2)
82  enddo
83  !$omp end do
84  !$omp end parallel
85 
86  initialized = .true.
87  hecmat%Iarray(98) = 0 ! symbolic setup done
88  hecmat%Iarray(97) = 0 ! numerical setup done
89 
90  end subroutine hecmw_precond_diag_22_setup
91 
92  subroutine hecmw_precond_diag_22_apply(WW)
93  implicit none
94  real(kind=kreal), intent(inout) :: ww(:)
95  integer(kind=kint) :: i
96  real(kind=kreal) :: x1, x2
97 
98  !C
99  !C== Block SCALING
100 
101  !$omp parallel default(none),private(i,X1,X2),shared(N,WW,ALU)
102  !$omp do
103  do i= 1, n
104  x1= ww(2*i-1)
105  x2= ww(2*i-0)
106  x2= x2 - alu(4*i-1)*x1
107  x2= alu(4*i )* x2
108  x1= alu(4*i-3)*( x1 - alu(4*i-2)*x2 )
109  ww(2*i-1)= x1
110  ww(2*i-0)= x2
111  enddo
112  !$omp end do
113  !$omp end parallel
114 
115  end subroutine hecmw_precond_diag_22_apply
116 
117  subroutine hecmw_precond_diag_22_clear()
118  implicit none
119  if (associated(alu)) deallocate(alu)
120  nullify(alu)
121  initialized = .false.
122  end subroutine hecmw_precond_diag_22_clear
123 
124 end module hecmw_precond_diag_22
hecmw_precond_diag_22::hecmw_precond_diag_22_setup
subroutine, public hecmw_precond_diag_22_setup(hecMAT)
Definition: hecmw_precond_DIAG_22.f90:28
hecmw_precond_diag_22
Definition: hecmw_precond_DIAG_22.f90:11
hecmw_precond_diag_22::hecmw_precond_diag_22_clear
subroutine, public hecmw_precond_diag_22_clear()
Definition: hecmw_precond_DIAG_22.f90:118
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_precond_diag_22::hecmw_precond_diag_22_apply
subroutine, public hecmw_precond_diag_22_apply(WW)
Definition: hecmw_precond_DIAG_22.f90:93
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_matrix_misc::hecmw_mat_get_sigma_diag
real(kind=kreal) function, public hecmw_mat_get_sigma_diag(hecMAT)
Definition: hecmw_matrix_misc.f90:694
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444