FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_precond_DIAG_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 
6 !C
7 !C***
8 !C*** module hecmw_precond_DIAG_11
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_11_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(1,1), pw(1)
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(np))
49  alu = 0.d0
50 
51  do ii= 1, n
52  alu(ii) = d(ii)
53  enddo
54 
55  !$omp parallel default(none),private(ii,ALUtmp,k,i,j,PW),shared(N,ALU,SIGMA_DIAG)
56  !$omp do
57  do ii= 1, n
58  alutmp(1,1)= alu(ii) * sigma_diag
59  alutmp(1,1)= 1.d0/alutmp(1,1)
60  alu(ii)= alutmp(1,1)
61  enddo
62  !$omp end do
63  !$omp end parallel
64 
65  initialized = .true.
66  hecmat%Iarray(98) = 0 ! symbolic setup done
67  hecmat%Iarray(97) = 0 ! numerical setup done
68 
69  end subroutine hecmw_precond_diag_11_setup
70 
71  subroutine hecmw_precond_diag_11_apply(WW)
72  implicit none
73  real(kind=kreal), intent(inout) :: ww(:)
74  integer(kind=kint) :: i
75 
76  !C
77  !C== Block SCALING
78 
79  !$omp parallel default(none),private(i),shared(N,WW,ALU)
80  !$omp do
81  do i= 1, n
82  ww(i)= alu(i)*ww(i)
83  enddo
84  !$omp end do
85  !$omp end parallel
86 
87  end subroutine hecmw_precond_diag_11_apply
88 
89  subroutine hecmw_precond_diag_11_clear()
90  implicit none
91  if (associated(alu)) deallocate(alu)
92  nullify(alu)
93  initialized = .false.
94  end subroutine hecmw_precond_diag_11_clear
95 
96 end module hecmw_precond_diag_11
hecmw_precond_diag_11::hecmw_precond_diag_11_apply
subroutine, public hecmw_precond_diag_11_apply(WW)
Definition: hecmw_precond_DIAG_11.f90:72
hecmw_precond_diag_11::hecmw_precond_diag_11_clear
subroutine, public hecmw_precond_diag_11_clear()
Definition: hecmw_precond_DIAG_11.f90:90
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
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_precond_diag_11
Definition: hecmw_precond_DIAG_11.f90:11
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_precond_diag_11::hecmw_precond_diag_11_setup
subroutine, public hecmw_precond_diag_11_setup(hecMAT)
Definition: hecmw_precond_DIAG_11.f90:28
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444