Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/fhutiter/src/huti_cg.F90
3206 views
1
module huti_cg
2
use huti_aux
3
implicit none
4
!
5
! *
6
! * Elmer, A Finite Element Software for Multiphysical Problems
7
! *
8
! * Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland
9
! *
10
! * This library is free software; you can redistribute it and/or
11
! * modify it under the terms of the GNU Lesser General Public
12
! * License as published by the Free Software Foundation; either
13
! * version 2.1 of the License, or (at your option) any later version.
14
! *
15
! * This library is distributed in the hope that it will be useful,
16
! * but WITHOUT ANY WARRANTY; without even the implied warranty of
17
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18
! * Lesser General Public License for more details.
19
! *
20
! * You should have received a copy of the GNU Lesser General Public
21
! * License along with this library (in file ../LGPL-2.1); if not, write
22
! * to the Free Software Foundation, Inc., 51 Franklin Street,
23
! * Fifth Floor, Boston, MA 02110-1301 USA
24
! *
25
! *****************************************************************************/
26
27
!
28
! Subroutines to implement Conjugate Gradient iterative method
29
!
30
31
#include "huti_fdefs.h"
32
33
!*************************************************************************
34
!*************************************************************************
35
!
36
! These subroutines are based on a book by Barret et al.:
37
! "Templates for the Solution of Linear Systems: Building Blocks for
38
! Iterative Methods", 1993.
39
!
40
! All matrix-vector operations are done externally, so we do not need
41
! to know about the matrix structure (sparse or dense). So has the
42
! memory allocation for the working arrays done also externally.
43
44
!*************************************************************************
45
! Work array is used in the following order:
46
! work(:,1) = z
47
! work(:,2) = p
48
! work(:,3) = q
49
! work(:,4) = r
50
!
51
! Definitions to make the code more understandable and to make it look
52
! like the pseudo code (these are common to all precisions)
53
!
54
55
#define X xvec
56
#define B rhsvec
57
58
#define Z work(:,1)
59
#define Z_ind 1
60
#define P work(:,2)
61
#define P_ind 2
62
#define Q work(:,3)
63
#define Q_ind 3
64
#define R work(:,4)
65
#define R_ind 4
66
!*************************************************************************
67
68
contains
69
70
71
!*************************************************************************
72
!*************************************************************************
73
! Single precision version
74
!*************************************************************************
75
!*************************************************************************
76
77
recursive subroutine huti_scgsolv ( ndim, wrkdim, xvec, rhsvec, &
78
ipar, dpar, work, matvecsubr, pcondlsubr, &
79
pcondrsubr, dotprodfun, normfun, stopcfun )
80
81
use huti_interfaces
82
implicit none
83
84
procedure( mv_iface_s ), pointer :: matvecsubr
85
procedure( pc_iface_s ), pointer :: pcondlsubr
86
procedure( pc_iface_s ), pointer :: pcondrsubr
87
procedure( dotp_iface_s ), pointer :: dotprodfun
88
procedure( norm_iface_s ), pointer :: normfun
89
procedure( stopc_iface_s ), pointer :: stopcfun
90
91
! Parameters
92
93
integer :: ndim, wrkdim
94
real, dimension(ndim) :: xvec, rhsvec
95
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
96
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
97
real, dimension(ndim,wrkdim) :: work
98
99
! Local variables
100
101
real :: alpha, beta, rho, oldrho
102
integer iter_count
103
real :: residual, rhsnorm, precrhsnorm
104
105
!
106
! End of variable declarations
107
!*********************************************************************
108
109
!*********************************************************************
110
! The actual CG begins here (look the pseudo code in the
111
! "Templates..."-book, page 15)
112
!
113
! First the initialization part
114
!
115
116
iter_count = 1
117
118
! Norms of right-hand side vector are used in convergence tests
119
120
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
121
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
122
rhsnorm = normfun( HUTI_NDIM, B, 1 )
123
end if
124
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
125
call pcondlsubr( P, B, ipar )
126
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
127
end if
128
129
! The following applies for all matrix operations in this solver
130
131
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
132
133
! Generate vector X if needed
134
135
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
136
call huti_srandvec ( X, ipar )
137
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
138
X = 1
139
end if
140
141
call matvecsubr( X, R, ipar )
142
143
R = B - R
144
145
!
146
! This is where the loop starts (that is we continue from here after
147
! the first iteration)
148
!
149
150
300 continue
151
152
call pcondlsubr( Q, R, ipar )
153
call pcondrsubr( Z, Q, ipar )
154
155
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
156
if ( rho .eq. 0 ) then
157
HUTI_INFO = HUTI_CG_RHO
158
go to 1000
159
end if
160
161
if ( iter_count .eq. 1 ) then
162
P = Z
163
else
164
beta = rho / oldrho
165
P = Z + beta * P
166
end if
167
168
call matvecsubr( P, Q, ipar )
169
170
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
171
172
X = X + alpha * P
173
R = R - alpha * Q
174
175
!
176
! Check the convergence against selected stopping criterion
177
!
178
179
select case (HUTI_STOPC)
180
case (HUTI_TRUERESIDUAL)
181
call matvecsubr( X, Z, ipar )
182
Z = Z - B
183
residual = normfun( HUTI_NDIM, Z, 1 )
184
case (HUTI_TRESID_SCALED_BYB)
185
call matvecsubr( X, Z, ipar )
186
Z = Z - B
187
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
188
case (HUTI_PSEUDORESIDUAL)
189
residual = normfun( HUTI_NDIM, R, 1 )
190
case (HUTI_PRESID_SCALED_BYB)
191
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
192
case (HUTI_PRESID_SCALED_BYPRECB)
193
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
194
case (HUTI_XDIFF_NORM)
195
Z = alpha * P
196
residual = normfun( HUTI_NDIM, Z, 1 )
197
case (HUTI_USUPPLIED_STOPC)
198
residual = stopcfun( X, B, R, ipar, dpar )
199
case default
200
call matvecsubr( X, Z, ipar )
201
Z = Z - B
202
residual = normfun( HUTI_NDIM, Z, 1 )
203
end select
204
205
!
206
! Print debugging output if desired
207
!
208
209
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
210
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
211
WRITE (*, '(I8, E11.4)') iter_count, residual
212
CALL FLUSH(6)
213
end if
214
end if
215
216
IF ( residual .LT. HUTI_TOLERANCE ) THEN
217
HUTI_INFO = HUTI_CONVERGENCE
218
go to 1000
219
end if
220
221
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
222
HUTI_INFO = HUTI_DIVERGENCE
223
GOTO 1000
224
END IF
225
226
oldrho = rho
227
228
!
229
! Return back to the iteration loop (without initialization)
230
!
231
232
iter_count = iter_count + 1
233
if ( iter_count .gt. HUTI_MAXIT ) then
234
HUTI_INFO = HUTI_MAXITER
235
go to 1000
236
end if
237
238
go to 300
239
240
!
241
! This is where we exit last time (after enough iterations or breakdown)
242
!
243
244
1000 continue
245
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
246
WRITE (*, '(I8, E11.4)') iter_count, residual
247
CALL FLUSH(6)
248
end if
249
250
HUTI_ITERS = iter_count
251
return
252
253
! End of execution
254
!*********************************************************************
255
256
end subroutine huti_scgsolv
257
258
!*************************************************************************
259
260
261
!*************************************************************************
262
!*************************************************************************
263
! Double precision version
264
!*************************************************************************
265
!*************************************************************************
266
267
recursive subroutine huti_dcgsolv ( ndim, wrkdim, xvec, rhsvec, &
268
ipar, dpar, work, matvecsubr, pcondlsubr, &
269
pcondrsubr, dotprodfun, normfun, stopcfun )
270
271
use huti_interfaces
272
implicit none
273
274
procedure( mv_iface_d ), pointer :: matvecsubr
275
procedure( pc_iface_d ), pointer :: pcondlsubr
276
procedure( pc_iface_d ), pointer :: pcondrsubr
277
procedure( dotp_iface_d ), pointer :: dotprodfun
278
procedure( norm_iface_d ), pointer :: normfun
279
procedure( stopc_iface_d ), pointer :: stopcfun
280
281
! Parameters
282
283
integer :: ndim, wrkdim
284
double precision, dimension(ndim) :: xvec, rhsvec
285
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
286
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
287
double precision, dimension(ndim,wrkdim) :: work
288
289
! Local variables
290
291
double precision :: alpha, beta, rho, oldrho
292
integer iter_count, i
293
double precision :: residual, rhsnorm, precrhsnorm
294
295
!
296
! End of variable declarations
297
!*********************************************************************
298
299
!*********************************************************************
300
! The actual CG begins here (look the pseudo code in the
301
! "Templates..."-book, page 15)
302
!
303
! First the initialization part
304
!
305
306
iter_count = 1
307
308
! Norms of right-hand side vector are used in convergence tests
309
310
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
311
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
312
rhsnorm = normfun( HUTI_NDIM, B, 1 )
313
end if
314
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
315
call pcondlsubr( P, B, ipar )
316
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
317
end if
318
319
! The following applies for all matrix operations in this solver
320
321
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
322
323
! Generate vector X if needed
324
325
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
326
call huti_drandvec ( X, ipar )
327
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
328
X = 1
329
end if
330
331
call matvecsubr( X, R, ipar )
332
333
#ifdef _OPENMP
334
!$OMP PARALLEL DO
335
do i=1,ndim
336
work(i,R_ind)=rhsvec(i)-work(i,R_ind)
337
end do
338
!$OMP END PARALLEL DO
339
#else
340
R = B - R
341
#endif
342
343
!
344
! This is where the loop starts (that is we continue from here after
345
! the first iteration)
346
!
347
348
300 continue
349
350
call pcondlsubr( Q, R, ipar )
351
call pcondrsubr( Z, Q, ipar )
352
353
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
354
if ( rho .eq. 0 ) then
355
HUTI_INFO = HUTI_CG_RHO
356
go to 1000
357
end if
358
359
if ( iter_count .eq. 1 ) then
360
#ifdef _OPENMP
361
!$OMP PARALLEL DO
362
DO i=1,ndim
363
work(i,P_ind)=work(i,Z_ind)
364
END DO
365
!$OMP END PARALLEL DO
366
#else
367
P = Z
368
#endif
369
else
370
beta = rho / oldrho
371
#ifdef _OPENMP
372
!$OMP PARALLEL DO
373
DO i=1,ndim
374
work(i,P_ind)=work(i,Z_ind) + beta * work(i,P_ind)
375
END DO
376
!$OMP END PARALLEL DO
377
#else
378
P = Z + beta * P
379
#endif
380
end if
381
382
call matvecsubr( P, Q, ipar )
383
384
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
385
386
#ifdef _OPENMP
387
!$OMP PARALLEL PRIVATE(i)
388
!$OMP DO
389
DO i=1,ndim
390
xvec(i) = xvec(i) + alpha * work(i,P_ind)
391
END DO
392
!$OMP END DO NOWAIT
393
!$OMP DO
394
DO i=1,ndim
395
work(i,R_ind) = work(i,R_ind) - alpha * work(i,Q_ind)
396
END DO
397
!$OMP END DO NOWAIT
398
!$OMP END PARALLEL
399
#else
400
X = X + alpha * P
401
R = R - alpha * Q
402
#endif
403
404
!
405
! Check the convergence against selected stopping criterion
406
!
407
408
select case (HUTI_STOPC)
409
case (HUTI_TRUERESIDUAL)
410
call matvecsubr( X, Z, ipar )
411
#ifdef _OPENMP
412
!$OMP PARALLEL DO
413
DO i=1,ndim
414
work(i,Z_ind) = work(i,Z_ind) - rhsvec(i)
415
END DO
416
!$OMP END PARALLEL DO
417
#else
418
Z = Z - B
419
#endif
420
residual = normfun( HUTI_NDIM, Z, 1 )
421
case (HUTI_TRESID_SCALED_BYB)
422
call matvecsubr( X, Z, ipar )
423
#ifdef _OPENMP
424
!$OMP PARALLEL DO
425
DO i=1,ndim
426
work(i,Z_ind) = work(i,Z_ind) - rhsvec(i)
427
END DO
428
!$OMP END PARALLEL DO
429
#else
430
Z = Z - B
431
#endif
432
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
433
case (HUTI_PSEUDORESIDUAL)
434
residual = normfun( HUTI_NDIM, R, 1 )
435
case (HUTI_PRESID_SCALED_BYB)
436
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
437
case (HUTI_PRESID_SCALED_BYPRECB)
438
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
439
case (HUTI_XDIFF_NORM)
440
#ifdef _OPENMP
441
!$OMP PARALLEL DO
442
DO i=1,ndim
443
work(i,Z_ind) = alpha * work(i,P_ind)
444
END DO
445
!$OMP END PARALLEL DO
446
#else
447
Z = alpha * P
448
#endif
449
residual = normfun( HUTI_NDIM, Z, 1 )
450
case (HUTI_USUPPLIED_STOPC)
451
residual = stopcfun( X, B, R, ipar, dpar )
452
case default
453
call matvecsubr( X, Z, ipar )
454
#ifdef _OPENMP
455
!$OMP PARALLEL DO
456
DO i=1,ndim
457
work(i,Z_ind) = work(i,Z_ind) - rhsvec(i)
458
END DO
459
!$OMP END PARALLEL DO
460
#else
461
Z = Z - B
462
#endif
463
residual = normfun( HUTI_NDIM, Z, 1 )
464
end select
465
466
!
467
! Print debugging output if desired
468
!
469
470
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
471
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
472
WRITE (*, '(I8, E11.4)') iter_count, residual
473
CALL FLUSH(6)
474
end if
475
end if
476
477
if ( residual .lt. HUTI_TOLERANCE ) then
478
HUTI_INFO = HUTI_CONVERGENCE
479
go to 1000
480
end if
481
482
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
483
HUTI_INFO = HUTI_DIVERGENCE
484
GOTO 1000
485
END IF
486
487
488
oldrho = rho
489
490
!
491
! Return back to the iteration loop (without initialization)
492
!
493
494
iter_count = iter_count + 1
495
if ( iter_count .gt. HUTI_MAXIT ) then
496
HUTI_INFO = HUTI_MAXITER
497
go to 1000
498
end if
499
500
go to 300
501
502
!
503
! This is where we exit last time (after enough iterations or breakdown)
504
!
505
506
1000 continue
507
IF ( HUTI_DBUGLVL .NE. HUTI_NO_DEBUG ) THEN
508
WRITE (*, '(I8, E11.4)') iter_count, residual
509
CALL FLUSH(6)
510
END IF
511
512
HUTI_ITERS = iter_count
513
return
514
515
! End of execution
516
!*********************************************************************
517
518
end subroutine huti_dcgsolv
519
520
!*************************************************************************
521
522
523
!*************************************************************************
524
!*************************************************************************
525
! Complex version
526
!*************************************************************************
527
!*************************************************************************
528
529
recursive subroutine huti_ccgsolv ( ndim, wrkdim, xvec, rhsvec, &
530
ipar, dpar, work, matvecsubr, pcondlsubr, &
531
pcondrsubr, dotprodfun, normfun, stopcfun )
532
533
use huti_interfaces
534
implicit none
535
536
procedure( mv_iface_c ), pointer :: matvecsubr
537
procedure( pc_iface_c ), pointer :: pcondlsubr
538
procedure( pc_iface_c ), pointer :: pcondrsubr
539
procedure( dotp_iface_c ), pointer :: dotprodfun
540
procedure( norm_iface_c ), pointer :: normfun
541
procedure( stopc_iface_c ), pointer :: stopcfun
542
543
! Parameters
544
545
integer :: ndim, wrkdim
546
complex, dimension(ndim) :: xvec, rhsvec
547
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
548
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
549
complex, dimension(ndim,wrkdim) :: work
550
551
! Local variables
552
553
complex :: alpha, beta, rho, oldrho
554
integer iter_count
555
real :: residual, rhsnorm, precrhsnorm
556
557
!
558
! End of variable declarations
559
!*********************************************************************
560
561
!*********************************************************************
562
! The actual CG begins here (look the pseudo code in the
563
! "Templates..."-book, page 15)
564
!
565
! First the initialization part
566
!
567
568
iter_count = 1
569
570
! Norms of right-hand side vector are used in convergence tests
571
572
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
573
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
574
rhsnorm = normfun( HUTI_NDIM, B, 1 )
575
end if
576
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
577
call pcondlsubr( P, B, ipar )
578
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
579
end if
580
581
! The following applies for all matrix operations in this solver
582
583
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
584
585
! Generate vector X if needed
586
587
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
588
call huti_crandvec ( X, ipar )
589
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
590
X = 1
591
end if
592
593
call matvecsubr( X, R, ipar )
594
595
R = B - R
596
597
!
598
! This is where the loop starts (that is we continue from here after
599
! the first iteration)
600
!
601
602
300 continue
603
604
call pcondlsubr( Q, R, ipar )
605
call pcondrsubr( Z, Q, ipar )
606
607
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
608
if ( rho .eq. 0 ) then
609
HUTI_INFO = HUTI_CG_RHO
610
go to 1000
611
end if
612
613
if ( iter_count .eq. 1 ) then
614
P = Z
615
else
616
beta = rho / oldrho
617
P = Z + beta * P
618
end if
619
620
call matvecsubr( P, Q, ipar )
621
622
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
623
624
X = X + alpha * P
625
R = R - alpha * Q
626
627
!
628
! Check the convergence against selected stopping criterion
629
!
630
631
select case (HUTI_STOPC)
632
case (HUTI_TRUERESIDUAL)
633
call matvecsubr( X, Z, ipar )
634
Z = Z - B
635
residual = normfun( HUTI_NDIM, Z, 1 )
636
case (HUTI_TRESID_SCALED_BYB)
637
call matvecsubr( X, Z, ipar )
638
Z = Z - B
639
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
640
case (HUTI_PSEUDORESIDUAL)
641
residual = normfun( HUTI_NDIM, R, 1 )
642
case (HUTI_PRESID_SCALED_BYB)
643
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
644
case (HUTI_PRESID_SCALED_BYPRECB)
645
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
646
case (HUTI_XDIFF_NORM)
647
Z = alpha * P
648
residual = normfun( HUTI_NDIM, Z, 1 )
649
case (HUTI_USUPPLIED_STOPC)
650
residual = stopcfun( X, B, R, ipar, dpar )
651
case default
652
call matvecsubr( X, Z, ipar )
653
Z = Z - B
654
residual = normfun( HUTI_NDIM, Z, 1 )
655
end select
656
657
!
658
! Print debugging output if desired
659
!
660
661
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
662
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
663
WRITE (*, '(I8, E11.4)') iter_count, residual
664
CALL FLUSH(6)
665
end if
666
end if
667
668
if ( residual .lt. HUTI_TOLERANCE ) then
669
HUTI_INFO = HUTI_CONVERGENCE
670
go to 1000
671
end if
672
673
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
674
HUTI_INFO = HUTI_DIVERGENCE
675
GOTO 1000
676
END IF
677
678
679
oldrho = rho
680
681
!
682
! Return back to the iteration loop (without initialization)
683
!
684
685
iter_count = iter_count + 1
686
if ( iter_count .gt. HUTI_MAXIT ) then
687
HUTI_INFO = HUTI_MAXITER
688
go to 1000
689
end if
690
691
go to 300
692
693
!
694
! This is where we exit last time (after enough iterations or breakdown)
695
!
696
697
1000 continue
698
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
699
WRITE (*, '(I8, E11.4)') iter_count, residual
700
call flush(6)
701
end if
702
703
HUTI_ITERS = iter_count
704
return
705
706
! End of execution
707
!*********************************************************************
708
709
end subroutine huti_ccgsolv
710
711
!*************************************************************************
712
713
714
!*************************************************************************
715
!*************************************************************************
716
! Double complex version
717
!*************************************************************************
718
!*************************************************************************
719
720
recursive subroutine huti_zcgsolv ( ndim, wrkdim, xvec, rhsvec, &
721
ipar, dpar, work, matvecsubr, pcondlsubr, &
722
pcondrsubr, dotprodfun, normfun, stopcfun )
723
724
use huti_interfaces
725
implicit none
726
727
procedure( mv_iface_z ), pointer :: matvecsubr
728
procedure( pc_iface_z ), pointer :: pcondlsubr
729
procedure( pc_iface_z ), pointer :: pcondrsubr
730
procedure( dotp_iface_z ), pointer :: dotprodfun
731
procedure( norm_iface_z ), pointer :: normfun
732
procedure( stopc_iface_z ), pointer :: stopcfun
733
734
! Parameters
735
736
integer :: ndim, wrkdim
737
double complex, dimension(ndim) :: xvec, rhsvec
738
integer, dimension(HUTI_IPAR_DFLTSIZE) :: ipar
739
double precision, dimension(HUTI_DPAR_DFLTSIZE) :: dpar
740
double complex, dimension(ndim,wrkdim) :: work
741
742
! Local variables
743
744
double complex :: alpha, beta, rho, oldrho
745
integer iter_count
746
double precision :: residual, rhsnorm, precrhsnorm
747
748
!
749
! End of variable declarations
750
!*********************************************************************
751
752
!*********************************************************************
753
! The actual CG begins here (look the pseudo code in the
754
! "Templates..."-book, page 15)
755
!
756
! First the initialization part
757
!
758
759
iter_count = 1
760
761
! Norms of right-hand side vector are used in convergence tests
762
763
if ( HUTI_STOPC .eq. HUTI_TRESID_SCALED_BYB .or. &
764
HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYB ) then
765
rhsnorm = normfun( HUTI_NDIM, B, 1 )
766
end if
767
if ( HUTI_STOPC .eq. HUTI_PRESID_SCALED_BYPRECB ) then
768
call pcondlsubr( P, B, ipar )
769
precrhsnorm = normfun( HUTI_NDIM, P, 1 )
770
end if
771
772
! The following applies for all matrix operations in this solver
773
774
HUTI_EXTOP_MATTYPE = HUTI_MAT_NOTTRPSED
775
776
! Generate vector X if needed
777
778
if ( HUTI_INITIALX .eq. HUTI_RANDOMX ) then
779
call huti_zrandvec ( X, ipar )
780
else if ( HUTI_INITIALX .ne. HUTI_USERSUPPLIEDX ) then
781
X = 1
782
end if
783
784
call matvecsubr( X, R, ipar )
785
786
R = B - R
787
788
!
789
! This is where the loop starts (that is we continue from here after
790
! the first iteration)
791
!
792
793
300 continue
794
795
call pcondlsubr( Q, R, ipar )
796
call pcondrsubr( Z, Q, ipar )
797
798
rho = dotprodfun( HUTI_NDIM, R, 1, Z, 1 )
799
if ( rho .eq. 0 ) then
800
HUTI_INFO = HUTI_CG_RHO
801
go to 1000
802
end if
803
804
if ( iter_count .eq. 1 ) then
805
P = Z
806
else
807
beta = rho / oldrho
808
P = Z + beta * P
809
end if
810
811
call matvecsubr( P, Q, ipar )
812
813
alpha = rho / dotprodfun( HUTI_NDIM, P, 1, Q, 1 )
814
815
X = X + alpha * P
816
R = R - alpha * Q
817
818
!
819
! Check the convergence against selected stopping criterion
820
!
821
822
select case (HUTI_STOPC)
823
case (HUTI_TRUERESIDUAL)
824
call matvecsubr( X, Z, ipar )
825
Z = Z - B
826
residual = normfun( HUTI_NDIM, Z, 1 )
827
case (HUTI_TRESID_SCALED_BYB)
828
call matvecsubr( X, Z, ipar )
829
Z = Z - B
830
residual = normfun( HUTI_NDIM, Z, 1 ) / rhsnorm
831
case (HUTI_PSEUDORESIDUAL)
832
residual = normfun( HUTI_NDIM, R, 1 )
833
case (HUTI_PRESID_SCALED_BYB)
834
residual = normfun( HUTI_NDIM, R, 1 ) / rhsnorm
835
case (HUTI_PRESID_SCALED_BYPRECB)
836
residual = normfun( HUTI_NDIM, R, 1 ) / precrhsnorm
837
case (HUTI_XDIFF_NORM)
838
Z = alpha * P
839
residual = normfun( HUTI_NDIM, Z, 1 )
840
case (HUTI_USUPPLIED_STOPC)
841
residual = stopcfun( X, B, R, ipar, dpar )
842
case default
843
call matvecsubr( X, Z, ipar )
844
Z = Z - B
845
residual = normfun( HUTI_NDIM, Z, 1 )
846
end select
847
848
!
849
! Print debugging output if desired
850
!
851
852
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
853
if ( mod(iter_count, HUTI_DBUGLVL) .eq. 0 ) then
854
WRITE (*, '(I8, E11.4)') iter_count, residual
855
call flush(6)
856
end if
857
end if
858
859
if ( residual .lt. HUTI_TOLERANCE ) then
860
HUTI_INFO = HUTI_CONVERGENCE
861
go to 1000
862
end if
863
864
IF( residual /= residual .OR. residual > HUTI_MAXTOLERANCE ) THEN
865
HUTI_INFO = HUTI_DIVERGENCE
866
GOTO 1000
867
END IF
868
869
oldrho = rho
870
871
!
872
! Return back to the iteration loop (without initialization)
873
!
874
875
iter_count = iter_count + 1
876
if ( iter_count .gt. HUTI_MAXIT ) then
877
HUTI_INFO = HUTI_MAXITER
878
go to 1000
879
end if
880
881
go to 300
882
883
!
884
! This is where we exit last time (after enough iterations or breakdown)
885
!
886
887
1000 continue
888
if ( HUTI_DBUGLVL .ne. HUTI_NO_DEBUG ) then
889
WRITE (*, '(I8, E11.4)') iter_count, residual
890
CALL FLUSH(6)
891
end if
892
893
HUTI_ITERS = iter_count
894
return
895
896
! End of execution
897
!*********************************************************************
898
899
end subroutine huti_zcgsolv
900
901
!*************************************************************************
902
903
904
end module huti_cg
905
906