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