SFEMaNS  version 5.3
Reference documentation for SFEMaNS
solver.f90
Go to the documentation of this file.
1 MODULE solve_petsc
2 USE my_util
4  INTEGER:: it_max
5  REAL(KIND=8) :: rel_tol, abs_tol
6  CHARACTER(LEN=20) :: solver, precond
7  LOGICAL :: verbose
8  END TYPE solver_param
9 
10 CONTAINS
11  SUBROUTINE init_solver(my_par,my_ksp,matrix,communicator,solver,precond, opt_re_init)
13 #include "petsc/finclude/petsc.h"
14  USE petsc
15  IMPLICIT NONE
16  LOGICAL, INTENT(IN), OPTIONAL :: opt_re_init
17  TYPE(solver_param) :: my_par
18  CHARACTER(*), OPTIONAL :: solver, precond
19  LOGICAL :: re_init
20  INTEGER :: deb, fin
21 
22  mat :: matrix
23  ksp :: my_ksp
24  pc :: prec
25  petscerrorcode :: ierr
26  mpi_comm :: communicator
27 
28  IF (.NOT.PRESENT(opt_re_init)) THEN
29  re_init=.false.
30  ELSE
31  re_init=opt_re_init
32  END IF
33 
34  IF (my_par%it_max.LE.0) THEN
35  my_par%it_max = 100
36  END IF
37  IF (my_par%rel_tol.LE.0.d0) THEN
38  my_par%rel_tol = 1.d-8
39  END IF
40  IF (my_par%abs_tol.LE.0.d0) THEN
41  my_par%abs_tol = 1.d-14
42  END IF
43 
44  IF (.NOT.re_init) CALL kspcreate(communicator,my_ksp,ierr)
45  !CALL KSPCreate(communicator,my_ksp,ierr)
46 
47  !CALL KSPSetOperators(my_ksp,matrix,matrix,DIFFERENT_NONZERO_PATTERN,ierr)
48  !CALL KSPSetOperators(my_ksp,matrix,matrix,SAME_NONZERO_PATTERN ,ierr) !Petsc 3.4.2
49  CALL kspsetoperators(my_ksp,matrix,matrix,ierr) !Petsc 3.7.2
50 
51  IF (PRESENT(solver)) THEN
52  deb = start_of_string(solver)
53  fin = last_of_string(solver)
54  IF (solver(deb:fin)=='BCGS') THEN
55  CALL kspsettype(my_ksp, kspbcgs, ierr)
56  ELSE IF (solver(deb:fin)=='GMRES') THEN
57  CALL kspsettype(my_ksp, kspgmres, ierr)
58  ELSE IF (solver(deb:fin)=='FGMRES') THEN
59  CALL kspsettype(my_ksp, kspfgmres, ierr)
60  ELSE IF (solver(deb:fin)=='PCR') THEN
61  CALL kspsettype(my_ksp, kspcr, ierr)
62  ELSE IF (solver(deb:fin)=='CHEBYCHEV') THEN
63  CALL kspsettype(my_ksp, kspchebyshev, ierr)
64  ELSE IF (solver(deb:fin)=='CG') THEN
65  CALL kspsettype(my_ksp, kspcg, ierr)
66  ELSE
67  CALL kspsettype(my_ksp, kspfgmres, ierr)
68  END IF
69  ELSE
70  CALL kspsettype(my_ksp, kspfgmres, ierr)
71  END IF
72  CALL kspsettolerances(my_ksp, my_par%rel_tol, my_par%abs_tol, &
73  petsc_default_real, my_par%it_max, ierr)
74  CALL kspgetpc(my_ksp, prec, ierr)
75  IF (PRESENT(precond)) THEN
76  deb = start_of_string(precond)
77  fin = last_of_string(precond)
78  IF (precond(deb:fin)=='JACOBI') THEN
79  CALL pcsettype(prec, pcbjacobi, ierr)
80  ELSE IF (precond(deb:fin)=='HYPRE') THEN
81  CALL pcsettype(prec, pchypre, ierr)
82  ELSE IF (precond(deb:fin)=='SSOR') THEN
83  CALL pcsettype(prec, pcsor, ierr)
84  ELSE IF (precond(deb:fin)=='MUMPS') THEN
85  CALL pcsettype(prec, pclu, ierr)
86  CALL kspsettype(my_ksp, ksppreonly, ierr)
87  !CALL PCFactorSetMatSolverPackage(prec, MATSOLVERMUMPS, ierr) !(JLG) Feb 20, 2019, petsc.3.8.4
88  CALL pcfactorsetmatsolvertype(prec, matsolvermumps, ierr) !
89  ELSE
90  CALL pcsettype(prec, pchypre, ierr)
91  END IF
92  ELSE
93  CALL pcsettype(prec, pchypre, ierr)
94  END IF
95  CALL kspsetfromoptions(my_ksp, ierr)
96  END SUBROUTINE init_solver
97 
98  SUBROUTINE solver(my_ksp,b,x,reinit,verbose)
99 #include "petsc/finclude/petsc.h"
100  use petsc
101  IMPLICIT NONE
102  LOGICAL, OPTIONAL :: reinit, verbose
103  INTEGER :: its
104  ksp :: my_ksp
105  petscerrorcode :: ierr
106  vec :: x, b
107  kspconvergedreason :: reason
108  IF (.NOT.PRESENT(reinit)) reinit=.true.
109  IF (.NOT.PRESENT(verbose)) verbose=.false.
110 
111  IF (reinit) CALL veczeroentries (x,ierr)
112  CALL kspsolve(my_ksp,b,x,ierr)
113  IF (verbose) THEN
114  CALL kspgetiterationnumber(my_ksp, its, ierr)
115  CALL kspgetconvergedreason(my_ksp, reason, ierr)
116  SELECT CASE(reason)
117  CASE(2)
118  WRITE(*,*) "KSP_CONVERGED_RTOL, Nb of iterations", its
119  CASE(3)
120  WRITE(*,*) "KSP_CONVERGED_ATOL, Nb of iterations", its
121  CASE(4)
122  WRITE(*,*) "Converged after one single iteration of the preconditioner is applied"
123  CASE(5,6,7,8)
124  WRITE(*,*) "Converge for strange reason:", reason
125  CASE(-2)
126  WRITE(*,*) "KSP_DIVERGED_NULL"
127  CASE(-3)
128  WRITE(*,*) "Not converged after it_max", its
129  CASE(-4)
130  WRITE(*,*) "Not converged: explosion"
131  CASE(-5,-6,-7)
132  WRITE(*,*) "Not converged for strange reasons", reason
133  CASE(-8)
134  WRITE(*,*) "Not converged: Indefinite preconditioner"
135  CASE(-9)
136  WRITE(*,*) "Not converged: NAN"
137  CASE(-10)
138  WRITE(*,*) "Not converged: Indefinite matrix"
139  CASE DEFAULT
140  WRITE(*,*) "Something strange happened", reason
141  END SELECT
142  END IF
143 
144  END SUBROUTINE solver
145 
146  SUBROUTINE create_local_petsc_matrix(communicator, LA, matrix, clean)
148 #include "petsc/finclude/petsc.h"
149  use petsc
150  IMPLICIT NONE
151  TYPE(petsc_csr_la) :: LA
152  LOGICAL, OPTIONAL :: clean
153  REAL(KIND=8),DIMENSION(:), POINTER :: aa
154  INTEGER :: nnzm1, dom_np
155  LOGICAL :: test_clean
156 !!$ INTEGER, DIMENSION(:), POINTER :: ia, ja
157  mpi_comm :: communicator
158  mat :: matrix
159  petscerrorcode :: ierr
160  !------------------------------------------------------------------------------
161  dom_np = SIZE(la%ia)-1
162  nnzm1=la%ia(dom_np)-la%ia(0)-1
163  ALLOCATE(aa(0:nnzm1))
164  aa =0.d0
165 
166 
167 !!$ ALLOCATE(ia(0:dom_np),ja(0:nnzm1))
168 !!$ ia = LA%ia
169 !!$ ja = LA%ja
170 !!$ CALL MatCreateMPIAIJWithArrays(communicator,dom_np,dom_np,PETSC_DECIDE, &
171 !!$ PETSC_DECIDE, ia, ja, aa, matrix, ierr)
172 !!$DEALLOCATE(ia,ja)
173 
174  CALL matcreatempiaijwitharrays(communicator,dom_np,dom_np,petsc_decide, &
175  petsc_decide, la%ia, la%ja, aa, matrix, ierr)
176 
177  DEALLOCATE(aa)
178  IF (PRESENT(clean)) THEN
179  test_clean=clean
180  ELSE
181  test_clean=.true.
182  END IF
183  IF (test_clean) THEN
184  IF (ASSOCIATED(la%ia)) DEALLOCATE(la%ia)
185  IF (ASSOCIATED(la%ja)) DEALLOCATE(la%ja)
186  END IF
187  END SUBROUTINE create_local_petsc_matrix
188 
189  SUBROUTINE create_local_petsc_matrix_a_detruire(communicator, aij, i_loc, matrix)
191 #include "petsc/finclude/petsc.h"
192  use petsc
193  IMPLICIT NONE
194  TYPE(aij_type), INTENT(IN) :: aij
195  INTEGER, DIMENSION(2) :: i_loc
196  INTEGER, DIMENSION(:), POINTER :: ia, ja
197  REAL(KIND=8),DIMENSION(:), POINTER :: aa
198  INTEGER :: nnzm1, dom_np, p, i, n
199 
200  mpi_comm :: communicator
201  mat :: matrix
202  petscerrorcode :: ierr
203  !------------------------------------------------------------------------------
204  dom_np = i_loc(2) - i_loc(1) + 1
205  nnzm1 = aij%ia(i_loc(2)+1)-aij%ia(i_loc(1))-1
206  ALLOCATE(ia(0:dom_np),ja(0:nnzm1))
207  ia(0)=0
208  DO i = 1, dom_np
209  n = i_loc(1) + i -1
210  ia(i) = aij%ia(n+1)-aij%ia(i_loc(1))
211  DO p=aij%ia(n), aij%ia(n+1)-1
212  ja(p-aij%ia(i_loc(1)))= aij%ja(p)-1
213  END DO
214  END DO
215  !------------------------------------------------------------------------------
216  ALLOCATE(aa(0:nnzm1))
217  aa =0
218  CALL matcreatempiaijwitharrays(communicator,dom_np,dom_np,petsc_decide, &
219  petsc_decide, ia, ja, aa, matrix, ierr)
220 
221  DEALLOCATE(ia,ja,aa)
223 
224  SUBROUTINE create_local_petsc_block_matrix(communicator, n_b, aij, i_loc, matrix)
226 #include "petsc/finclude/petsc.h"
227  use petsc
228  IMPLICIT NONE
229  TYPE(aij_type), INTENT(IN) :: aij
230  INTEGER, DIMENSION(2) :: i_loc
231  INTEGER :: n_b
232  INTEGER, DIMENSION(:), POINTER :: ia, ja
233  REAL(KIND=8),DIMENSION(:), POINTER :: aa
234  INTEGER :: nnzm1, dom_np, p, i, n, ib, k
235 
236  mpi_comm :: communicator
237  mat :: matrix
238  petscerrorcode :: ierr
239 
240  dom_np = i_loc(2) - i_loc(1) + 1
241  nnzm1 = n_b*(aij%ia(i_loc(2)+1)-aij%ia(i_loc(1))-1)
242  ALLOCATE(ia(0:n_b*dom_np),ja(0:nnzm1))
243  ia(0)=0
244  DO k = 1, n_b
245  DO i = 1, dom_np
246  ib = i + (k-1)*dom_np
247  n = i_loc(1) + i - 1
248  ia(i) = n_b*(aij%ia(n+1)-aij%ia(i_loc(1)))
249  DO p=aij%ia(n), aij%ia(n+1)-1
250  ja(p-aij%ia(i_loc(1)))= aij%ja(p)-1
251  END DO
252  END DO
253  END DO
254  !------------------------------------------------------------------------------
255  ALLOCATE(aa(0:nnzm1))
256  aa =0
257  CALL matcreatempiaijwitharrays(communicator,dom_np,dom_np,petsc_decide, &
258  petsc_decide, ia, ja, aa, matrix, ierr)
259 
260  DEALLOCATE(ia,ja,aa)
261  END SUBROUTINE create_local_petsc_block_matrix
262 
263 
264 END MODULE solve_petsc
subroutine solver(my_ksp, b, x, reinit, verbose)
Definition: solver.f90:99
integer function, public last_of_string(string)
subroutine create_local_petsc_matrix(communicator, LA, matrix, clean)
Definition: solver.f90:147
integer function, public start_of_string(string)
subroutine init_solver(my_par, my_ksp, matrix, communicator, solver, precond, opt_re_init)
Definition: solver.f90:12
subroutine create_local_petsc_block_matrix(communicator, n_b, aij, i_loc, matrix)
Definition: solver.f90:225
subroutine create_local_petsc_matrix_a_detruire(communicator, aij, i_loc, matrix)
Definition: solver.f90:190