FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_precond_DIAG_33.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_33
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_33_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(3,3), pw(3)
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(9*np))
49  alu = 0.d0
50 
51  do ii= 1, n
52  alu(9*ii-8) = d(9*ii-8)
53  alu(9*ii-7) = d(9*ii-7)
54  alu(9*ii-6) = d(9*ii-6)
55  alu(9*ii-5) = d(9*ii-5)
56  alu(9*ii-4) = d(9*ii-4)
57  alu(9*ii-3) = d(9*ii-3)
58  alu(9*ii-2) = d(9*ii-2)
59  alu(9*ii-1) = d(9*ii-1)
60  alu(9*ii ) = d(9*ii )
61  enddo
62 
63  !$omp parallel default(none),private(ii,ALUtmp,k,i,j,PW),shared(N,ALU,SIGMA_DIAG)
64  !$omp do
65  do ii= 1, n
66  alutmp(1,1)= alu(9*ii-8) * sigma_diag
67  alutmp(1,2)= alu(9*ii-7)
68  alutmp(1,3)= alu(9*ii-6)
69  alutmp(2,1)= alu(9*ii-5)
70  alutmp(2,2)= alu(9*ii-4) * sigma_diag
71  alutmp(2,3)= alu(9*ii-3)
72  alutmp(3,1)= alu(9*ii-2)
73  alutmp(3,2)= alu(9*ii-1)
74  alutmp(3,3)= alu(9*ii ) * sigma_diag
75  do k= 1, 3
76  alutmp(k,k)= 1.d0/alutmp(k,k)
77  do i= k+1, 3
78  alutmp(i,k)= alutmp(i,k) * alutmp(k,k)
79  do j= k+1, 3
80  pw(j)= alutmp(i,j) - alutmp(i,k)*alutmp(k,j)
81  enddo
82  do j= k+1, 3
83  alutmp(i,j)= pw(j)
84  enddo
85  enddo
86  enddo
87  alu(9*ii-8)= alutmp(1,1)
88  alu(9*ii-7)= alutmp(1,2)
89  alu(9*ii-6)= alutmp(1,3)
90  alu(9*ii-5)= alutmp(2,1)
91  alu(9*ii-4)= alutmp(2,2)
92  alu(9*ii-3)= alutmp(2,3)
93  alu(9*ii-2)= alutmp(3,1)
94  alu(9*ii-1)= alutmp(3,2)
95  alu(9*ii )= alutmp(3,3)
96  enddo
97  !$omp end do
98  !$omp end parallel
99 
100  initialized = .true.
101  hecmat%Iarray(98) = 0 ! symbolic setup done
102  hecmat%Iarray(97) = 0 ! numerical setup done
103 
104  end subroutine hecmw_precond_diag_33_setup
105 
106  subroutine hecmw_precond_diag_33_apply(WW)
107  implicit none
108  real(kind=kreal), intent(inout) :: ww(:)
109  integer(kind=kint) :: i
110  real(kind=kreal) :: x1, x2, x3
111 
112  !C
113  !C== Block SCALING
114 
115  !$omp parallel default(none),private(i,X1,X2,X3),shared(N,WW,ALU)
116  !$omp do
117  do i= 1, n
118  x1= ww(3*i-2)
119  x2= ww(3*i-1)
120  x3= ww(3*i )
121  x2= x2 - alu(9*i-5)*x1
122  x3= x3 - alu(9*i-2)*x1 - alu(9*i-1)*x2
123  x3= alu(9*i )* x3
124  x2= alu(9*i-4)*( x2 - alu(9*i-3)*x3 )
125  x1= alu(9*i-8)*( x1 - alu(9*i-6)*x3 - alu(9*i-7)*x2)
126  ww(3*i-2)= x1
127  ww(3*i-1)= x2
128  ww(3*i )= x3
129  enddo
130  !$omp end do
131  !$omp end parallel
132 
133  end subroutine hecmw_precond_diag_33_apply
134 
135  subroutine hecmw_precond_diag_33_clear()
136  implicit none
137  if (associated(alu)) deallocate(alu)
138  nullify(alu)
139  initialized = .false.
140  end subroutine hecmw_precond_diag_33_clear
141 
142 end module hecmw_precond_diag_33
hecmw_precond_diag_33::hecmw_precond_diag_33_apply
subroutine, public hecmw_precond_diag_33_apply(WW)
Definition: hecmw_precond_DIAG_33.f90:107
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_precond_diag_33::hecmw_precond_diag_33_clear
subroutine, public hecmw_precond_diag_33_clear()
Definition: hecmw_precond_DIAG_33.f90:136
hecmw_precond_diag_33
Definition: hecmw_precond_DIAG_33.f90:11
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_precond_diag_33::hecmw_precond_diag_33_setup
subroutine, public hecmw_precond_diag_33_setup(hecMAT)
Definition: hecmw_precond_DIAG_33.f90:28
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444