forked from ESCOMP/CLUBB_CESM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
gmres_wrap.F90
383 lines (313 loc) · 14.6 KB
/
gmres_wrap.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
!----------------------------------------------------------------------------
! $Id$
!==============================================================================
module gmres_wrap
#ifdef MKL
! Description:
! This module wraps the MKL version of GMRES, an iterative solver. Note that
! this will only work for the MKL-specific version of GMRES--any other GMRES
! implementations will require retooling of this code!
!
! The primary subroutine, gmres_solve utilizes GMRES to solve a given matrix.
!
! There is also a gmres_init, which initializes some of the internal data
! used for the wrapper.
!
! This wrapper automatically keeps prior solutions to use the previous data
! to speed up the solves. For the purposes of allowing this solver to be used
! with more than one matrix type, the wrapper has a "solve index" variable.
! Pass in the proper solve index variable to associate your solve with
! previous solves of the same matrix.
use gmres_cache, only: &
maximum_gmres_idx ! Variable
implicit none
public :: gmres_solve, gmres_init
private ! Default scope
contains
subroutine gmres_init(max_numeqns, max_elements) ! Intent(in)
! Description:
! Initialization subroutine for the GMRES iterative matrix equation solver
!
! This subroutine initializes the previous memory handles for the GMRES
! routines, for the purpose of speeding up calculations.
! These handles are initialized to a size specified by the number of
! equations specified in this subroutine.
!
! WARNING: Once initialized, only use the specified gmres_idx for that
! particular matrix! Failure to do so could result in greatly decreased
! performance, incorrect solutions, or both!
!
! Once this is called, the proper prev_soln_<matrix> and prev_lu_<matrix>
! handles in the gmres_cache module can be used, and will need to be passed
! in to gmres_solve for that matrix.
!
! References:
! None
use gmres_cache, only: &
gmres_cache_matrix_init ! Subroutines
implicit none
! Input Variables
integer, intent(in) :: &
max_numeqns, & ! Maximum number of equations for a matrix that will be
! solved with GMRES
max_elements ! Maximum number of non-zero elements for a matrix that
! will be solved with GMRES
call gmres_cache_matrix_init( max_numeqns, max_elements, maximum_gmres_idx )
end subroutine gmres_init
subroutine gmres_solve( elements, numeqns, & !Intent(in)
csr_a, csr_ia, csr_ja, tempsize, & !Intent(in)
prev_soln, prev_lu, rhs, temp, & !Intent(in/out)
solution ) !Intent(out)
! Description:
! Solves a matrix equation using GMRES. On the first timestep and every
! fifth timestep afterward, a preconditioner is computed for the matrix
! and stored. In addition, on the first timestep the matrix is solved using
! LAPACK, which is used as the estimate for GMRES for the first timestep.
! After this, the previous solution found is used as the estimate.
!
! To use the proper cached preconditioner and solution, make sure you pass
! the proper gmres_idx corresponding to the matrix you're solving--using a
! value different than what has been used in the past will cause, at best,
! a slower solve, and at worst, an incorrect one.
!
! References:
! None
use clubb_precision, only: &
dp, & ! double precision
core_rknd
use error_code, only: &
err_code ! Error indicator
implicit none
include "mkl_rci.fi"
! Input variables
integer, intent(in) :: &
elements, & ! Number of elements in the csr_a/csr_ja arrays
numeqns ! Number of equations in the matrix
real( kind = core_rknd ), dimension(elements), intent(in) :: &
csr_a ! A-array description of the matrix in CSR format. This
! will be converted to double precision for the purposes
! of running GMRES.
integer, dimension(numeqns + 1), intent(in) :: &
csr_ia ! IA-array portion of the matrix description in CSR format.
! This describes the indices of the JA-array that start
! new rows. For more details, check the documentation in
! the csr_matrix_module.
integer, dimension(elements), intent(in) :: &
csr_ja ! JA-array portion of the matrix description in CSR format.
! This describes which columns of a are nonzero. For more
! details, check the documentation in the csr_matrix_module.
integer, intent(in) :: &
tempsize ! Denotes the size of the temporary array used for GMRES
! calculations.
! Input/Output variables
real( kind = core_rknd ), dimension(numeqns), intent(inout) :: &
rhs ! Right-hand-side vectors to solve the equation for.
real( kind = dp ), dimension(numeqns), intent(inout) :: &
prev_soln ! Previous solution cache vector for the matrix to be solved
! for--pass the proper handle from the gmres_cache module
real( kind = dp ), dimension(elements), intent(inout) :: &
prev_lu ! Previous LU-decomposition a-array for the matrix to be
! solved for--pass the proper handle from the gmres_cache
! module
real( kind = dp ), dimension(tempsize), intent(inout) :: &
temp ! Temporary array that stores working values while the GMRES
! solver iterates
! Output variables
real( kind = core_rknd ), dimension(numeqns), intent(out) :: &
solution ! Solution vector, output of solver routine
! Local variables
logical :: l_gmres_run ! Variable denoting if we need to loop and run
! a GMRES iteration again.
integer :: &
rci_req, & ! RCI_Request for GMRES--allows us to take action based
! on what the iterative solver requests to be done.
iters ! Total number of iterations GMRES has run.
integer, dimension(128) :: &
ipar ! Parameter array for the GMRES iterative solver
real( kind = dp ), dimension(128) :: &
dpar ! Parameter array for the GMRES iterative solver
! The following local variables are double-precision so we can use GMRES
! as there is only double-precision support for GMRES.
! We will need to convert our single-precision numbers to double precision
! for the duration of the calculations.
real( kind = dp ), dimension(elements) :: &
csr_dbl_a ! Double-precision version of the CSR-format A array
real( kind = dp ), dimension(numeqns) :: &
dbl_rhs, & ! Double-precision version of the rhs vector
dbl_soln, & ! Double-precision version of the solution vector
tempvec ! Temporary vector for applying inverse LU-decomp matrix
!tmp_rhs
! Variables used to solve the preconditioner the first time with PARDISO.
!integer, parameter :: &
!pardiso_size_arrays = 64, &
!real_nonsymm = 11
!integer(kind=8), dimension(pardiso_size_arrays) :: &
! pt ! PARDISO internal pointer array
!integer(kind=4), dimension(pardiso_size_arrays) :: &
! iparm
!integer(kind=4), dimension(numeqns) :: &
! perm
! integer :: i, j
! We want to be running, initially.
l_gmres_run = .true.
! Convert our A array and rhs vector to double precision...
csr_dbl_a = real(csr_a, kind=dp)
dbl_rhs = real(rhs, kind=dp)
! DEBUG: Set our a_array so it represents the identity matrix, and
! set the RHS so we can get a meaningful answer.
! csr_dbl_a = 1_dp
! csr_dbl_a(1) = 1D1
! csr_dbl_a(5) = 1D1
! csr_dbl_a(elements) = 1D1
! csr_dbl_a(elements - 4) = 1D1
! do i=10,elements - 9,5
! csr_dbl_a(i) = 1D1
! end do
! do i=1,numeqns,1
! dbl_rhs(i) = i * 1_dp
! end do
! dbl_rhs = 9D3
! dbl_rhs = 1D1
! DEBUG: Make sure our a_array isn't wrong
! do i=1,elements,1
! print *, "csr_dbl_a idx",i,"=",csr_dbl_a(i)
! end do
! Figure out the default value for ipar(15) and put it in our ipar_15 int.
!ip_15 = min(150, numeqns)
! Figure out the size of the temp array.
!tempsize = ((((2*numeqns + 1)*numeqns)+(numeqns*(numeqns+9))/2) + 1)
! This ugly equation was lifted from the Intel documentation of dfgmres:
! http://www.intel.com/software/products/mkl/docs/webhelp/ssr/functn_rci_dfgmres.html
! All of the ipar(15)s have been replaced with "numeqns", as the code
! examples seemed to use N (numeqns) in place of ipar(15).
! Allocate the temp array.
!allocate(temp(1:tempsize))
! Generate our preconditioner matrix with the ILU0 subroutine.
call dcsrilu0( numeqns, csr_dbl_a, csr_ia, csr_ja, &
prev_lu, ipar, dpar )
! On the first timestep we need to solve our preconditioner to give us
! our first solution estimate. After this, the previous solution will
! suffice as an estimate.
! if (iteration_num(gmres_idx) == 0) then
!solve with precond_a, csr_ia, csr_ja.
!One thing to test, too: try just setting the solution vector to 1
! for the first timestep and see if it's not too unreasonably slow?
! call pardisoinit( pt, real_nonsymm, iparm )
#ifdef _OPENMP
! iparm(3) = omp_get_max_threads()
#else
! iparm(3) = 1
#endif
! call pardiso( pt, 1, 1, real_nonsymm, 13, numeqns, & !Intent(in)
! prev_lu, csr_ia, csr_ja, perm, 1, iparm, 0, & !Intent(in)
! dbl_rhs, & !Intent(inout)
! prev_soln ) !Intent(out)
! end if !iteration_num == 1
!DEBUG: Set apporximate solution vector to 0.9 (?) for now
!prev_soln(:) = 0.9_dp
!do i=1,numeqns,1
! print *, "Current approximate solution idx",i,"=",prev_soln(i)
!end do
! Initialize our solution vector to the previous solution passed in
dbl_soln = prev_soln
! Set up the GMRES solver.
call dfgmres_init( numeqns, dbl_soln, dbl_rhs, &
rci_req, ipar, dpar, temp )
! Set the parameters that tell GMRES to handle stopping tests
ipar(9) = 1
ipar(10) = 0
ipar(12) = 1
! Set the parameter that tells GMRES to use a preconditioner
ipar(11) = 1
! Check our GMRES settings.
call dfgmres_check( numeqns, dbl_soln, dbl_rhs, &
rci_req, ipar, dpar, temp )
! Start the GMRES solver. We set up a while loop which will be broken when
! the GMRES solver indicates that a solution has been found.
do while(l_gmres_run)
!print *, "********************************************************"
!print *, "BEGINNING ANOTHER ITERATION..."
!print *, "========================================================"
! Run a GMRES iteration.
call dfgmres( numeqns, dbl_soln, dbl_rhs, &
rci_req, ipar, dpar, temp )
select case(rci_req)
case (0)
l_gmres_run = .false.
case (1)
! Multiply our left-hand side by the vector placed in the temp array,
! at ipar(22), and place the result in the temp array at ipar(23).
! Display temp(ipar(22))
! print *, "------------------------------------------------"
! print *, "RCI_REQ=1: MULTIPLY VECTOR BY A MATRIX"
! do i=1,numeqns,1
! print *, "Tempvec before, idx",i,"=",temp(ipar(22)+i-1)
! end do
call mkl_dcsrgemv( 'N', numeqns, csr_dbl_a, csr_ia, csr_ja, &
temp(ipar(22)), temp(ipar(23)) ) ! Known magic number
! do i=1,numeqns,1
! print *, "Tempvec after, idx",i,"=",temp(ipar(23)+i-1)
! end do
! print *, "------------------------------------------------"
case (2)
! Ignore this for now, see if GMRES ever escapes.
case (3)
! Apply the inverse of the preconditioner to the vector placed in the
! temp array at ipar(22), and place the result in the temp array at
! ipar(23).
!print *, "------------------------------------------------"
!print *, "RCI_REQ=3: APPLY PRECONDITION TO VECTOR"
!do i=1,numeqns,1
! print *, "Tempvec before, idx",i,"=",temp(ipar(22)+i-1)
!end do
call mkl_dcsrtrsv( 'L', 'N', 'U', numeqns, &
prev_lu, csr_ia, csr_ja, &
temp(ipar(22)), tempvec ) ! Known magic number
call mkl_dcsrtrsv( 'U', 'N', 'N', numeqns, &
prev_lu, csr_ia, csr_ja, &
tempvec, temp(ipar(23)) ) ! Known magic number
!do i=1,numeqns,1
! print *, "Tempvec after, idx",i,"=",temp(ipar(23)+i-1)
!end do
!print *, "------------------------------------------------"
case (4)
! if (dpar(7) < GMRES_TOL) then
! l_gmres_run = .false.
! else
! ! Keep running, we aren't there yet.
! l_gmres_run = .true.
! end if
case default
! We got a response we weren't expecting. This is probably bad.
! (Then again, maybe it's just not something we accounted for?)
! Regardless, let's set an error code and break out of here.
write(fstderr,*) "Unknown rci_request returned from GMRES:", rci_req
err_code = clubb_fatal_error
return
end select
! Report current iteration
! call dfgmres_get( numeqns, dbl_soln, dbl_rhs, rci_req, &
! ipar, dpar, temp, iters )
! print *, "========================================================"
! print *, "END OF LOOP: REPORTING INFORMATION"
! print *, "Current number of GMRES iterations: ", iters
! do i=1,numeqns,1
! print *, "double value of soln so far, idx",i,"=",dbl_soln(i)
! end do
! print *, "========================================================"
! print *, "********************************************************"
end do
! Get the answer, convert it to single-precision
call dfgmres_get( numeqns, dbl_soln, dbl_rhs, rci_req, &
ipar, dpar, temp, iters )
!print *, "Total iterations for GMRES:",iters
!do i=1,numeqns,1
! print *, "double value of soln, idx",i,"=",dbl_soln(i)
!end do
! Store our solution as the previous solution for use in the next
! simulation timestep.
prev_soln = dbl_soln
solution = real(dbl_soln)
end subroutine gmres_solve
#endif /* MKL */
end module gmres_wrap