Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/fem/src/IterSolve.F90
3203 views
1
!/*****************************************************************************/
2
! *
3
! * Elmer, A Finite Element Software for Multiphysical Problems
4
! *
5
! * Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland
6
! *
7
! * This library is free software; you can redistribute it and/or
8
! * modify it under the terms of the GNU Lesser General Public
9
! * License as published by the Free Software Foundation; either
10
! * version 2.1 of the License, or (at your option) any later version.
11
! *
12
! * This library is distributed in the hope that it will be useful,
13
! * but WITHOUT ANY WARRANTY; without even the implied warranty of
14
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
15
! * Lesser General Public License for more details.
16
! *
17
! * You should have received a copy of the GNU Lesser General Public
18
! * License along with this library (in file ../LGPL-2.1); if not, write
19
! * to the Free Software Foundation, Inc., 51 Franklin Street,
20
! * Fifth Floor, Boston, MA 02110-1301 USA
21
! *
22
! *****************************************************************************/
23
!
24
! ******************************************************************************
25
! *
26
! * Authors: Juha Ruokolainen
27
! * Email: [email protected]
28
! * Web: http://www.csc.fi/elmer
29
! * Address: CSC - IT Center for Science Ltd.
30
! * Keilaranta 14
31
! * 02101 Espoo, Finland
32
! *
33
! * Original Date: 08 Jun 1997
34
! *
35
! ****************************************************************************/
36
37
#include "huti_fdefs.h"
38
39
!> \ingroup ElmerLib
40
!> \{
41
42
43
44
45
!------------------------------------------------------------------------------
46
!> Module containing control of the iterative solvers for linear systems
47
!> that come with the Elmer suite.
48
!------------------------------------------------------------------------------
49
MODULE IterSolve
50
51
USE Lists
52
USE BandMatrix
53
USE IterativeMethods
54
USE huti_sfe
55
56
IMPLICIT NONE
57
58
!/*
59
! * Iterative method selection
60
! */
61
INTEGER, PARAMETER, PRIVATE :: ITER_BiCGStab = 320
62
INTEGER, PARAMETER, PRIVATE :: ITER_TFQMR = 330
63
INTEGER, PARAMETER, PRIVATE :: ITER_CG = 340
64
INTEGER, PARAMETER, PRIVATE :: ITER_CGS = 350
65
INTEGER, PARAMETER, PRIVATE :: ITER_GMRES = 360
66
INTEGER, PARAMETER, PRIVATE :: ITER_BiCGStab2 = 370
67
INTEGER, PARAMETER, PRIVATE :: ITER_SGS = 380
68
INTEGER, PARAMETER, PRIVATE :: ITER_JACOBI = 390
69
INTEGER, PARAMETER, PRIVATE :: ITER_RICHARDSON = 391
70
INTEGER, PARAMETER, PRIVATE :: ITER_BICGSTABL = 400
71
INTEGER, PARAMETER, PRIVATE :: ITER_GCR = 410
72
INTEGER, PARAMETER, PRIVATE :: ITER_IDRS = 420
73
74
!/*
75
! * Preconditioning type code
76
! */
77
INTEGER, PARAMETER, PRIVATE :: PRECOND_NONE = 500
78
INTEGER, PARAMETER, PRIVATE :: PRECOND_DIAGONAL = 510
79
INTEGER, PARAMETER, PRIVATE :: PRECOND_ILUn = 520
80
INTEGER, PARAMETER, PRIVATE :: PRECOND_ILUT = 530
81
INTEGER, PARAMETER, PRIVATE :: PRECOND_MG = 540
82
INTEGER, PARAMETER, PRIVATE :: PRECOND_BILUn = 550
83
INTEGER, PARAMETER, PRIVATE :: PRECOND_Vanka = 560
84
INTEGER, PARAMETER, PRIVATE :: PRECOND_Circuit = 570
85
INTEGER, PARAMETER, PRIVATE :: PRECOND_Slave = 580
86
87
INTEGER, PARAMETER :: stack_max=64
88
INTEGER :: stack_pos=0
89
LOGICAL :: FirstCall(stack_max)
90
91
REAL(KIND=dp), POINTER :: fm_Diag(:), fm_G(:,:)
92
93
CONTAINS
94
95
96
#ifndef HUTI_MAXTOLERANCE
97
#define HUTI_MAXTOLERANCE dpar(2)
98
#endif
99
#ifndef HUTI_SGSPARAM
100
#define HUTI_SGSPARAM dpar(3)
101
#endif
102
#ifndef HUTI_PSEUDOCOMPLEX
103
#define HUTI_PSEUDOCOMPLEX ipar(7)
104
#endif
105
#ifndef HUTI_BICGSTABL_L
106
#define HUTI_BICGSTABL_L ipar(16)
107
#endif
108
#ifndef HUTI_DIVERGENCE
109
#define HUTI_DIVERGENCE 3
110
#endif
111
#ifndef HUTI_GCR_RESTART
112
#define HUTI_GCR_RESTART ipar(17)
113
#endif
114
#ifndef HUTI_IDRS_S
115
#define HUTI_IDRS_S ipar(18)
116
#endif
117
118
!------------------------------------------------------------------------------
119
!> Dummy preconditioner, if linear system scaling is active this corresponds
120
!> to diagonal preconditioning.
121
!------------------------------------------------------------------------------
122
SUBROUTINE pcond_dummy(u,v,ipar )
123
!------------------------------------------------------------------------------
124
INTEGER :: ipar(*)
125
REAL(KIND=dp) :: u(HUTI_NDIM), v(HUTI_NDIM)
126
INTEGER :: i
127
!------------------------------------------------------------------------------
128
!$OMP PARALLEL DO
129
DO i=1,HUTI_NDIM
130
u(i) = v(i)
131
END DO
132
!$OMP END PARALLEL DO
133
!------------------------------------------------------------------------------
134
END SUBROUTINE pcond_dummy
135
!------------------------------------------------------------------------------
136
137
!------------------------------------------------------------------------------
138
!> Complex dummy preconditioner, if linear system scaling is active this corresponds
139
!> to diagonal preconditioning.
140
!------------------------------------------------------------------------------
141
SUBROUTINE pcond_dummy_cmplx(u,v,ipar )
142
!------------------------------------------------------------------------------
143
INTEGER :: ipar(*)
144
COMPLEX(KIND=dp) :: u(HUTI_NDIM), v(HUTI_NDIM)
145
!------------------------------------------------------------------------------
146
u = v
147
!------------------------------------------------------------------------------
148
END SUBROUTINE pcond_dummy_cmplx
149
!------------------------------------------------------------------------------
150
151
!------------------------------------------------------------------------------
152
SUBROUTINE fm_DiagPrec( u,v,ipar )
153
!------------------------------------------------------------------------------
154
IMPLICIT NONE
155
156
REAL(KIND=dp) :: u(*),v(*)
157
INTEGER :: ipar(*)
158
159
INTEGER :: n
160
161
n = HUTI_NDIM
162
u(1:n) = v(1:n)*fm_diag(1:n)
163
!------------------------------------------------------------------------------
164
END SUBROUTINE fm_DiagPrec
165
!------------------------------------------------------------------------------
166
167
168
!------------------------------------------------------------------------------
169
SUBROUTINE fm_MatVec( u,v,ipar )
170
!------------------------------------------------------------------------------
171
IMPLICIT NONE
172
173
INTEGER :: ipar(*)
174
REAL(KIND=dp) :: u(*),v(*), ct, rsum, cumt=0, s
175
176
INTEGER :: i,j,n
177
178
n = HUTI_NDIM
179
#if 1
180
! CALL DSYMV('U',n,1.0_dp,Jacobian,n,u,1,0.0_dp,v,1)
181
CALL DGEMV('N',n,n,1.0_dp,fm_G,n,u,1,0.0_dp,v,1)
182
#else
183
v(1:n) = 0
184
!$omp parallel do private(i,j,s) shared(n,u,v,fM_G)
185
DO i=1,n
186
s = 0._dp
187
DO j=1,n
188
s = s + fm_G(i,j) * u(j)
189
END DO
190
v(i) = s
191
END DO
192
!$omp end parallel do
193
#endif
194
!------------------------------------------------------------------------------
195
END SUBROUTINE fm_MatVec
196
!------------------------------------------------------------------------------
197
198
199
!> \}
200
!> \}
201
202
!------------------------------------------------------------------------------
203
!> The routine that decides which linear system solver to call, and calls it.
204
!> There are two main sources of iterations within Elmer.
205
!> 1) The old HUTiter library that includes the most classic iterative Krylov
206
!> methods.
207
!> 2) The internal MODULE IterativeMethods that includes some classic iterative
208
!> methods and also some more recent Krylov methods.
209
!------------------------------------------------------------------------------
210
RECURSIVE SUBROUTINE IterSolver( A,x,b,Solver,ndim,DotF, &
211
NormF,MatvecF,PrecF,StopcF )
212
!------------------------------------------------------------------------------
213
USE huti_sfe
214
USE ListMatrix
215
USE SParIterGlobals
216
IMPLICIT NONE
217
218
!------------------------------------------------------------------------------
219
TYPE(Solver_t) :: Solver
220
REAL(KIND=dp), DIMENSION(:), TARGET CONTIG :: x,b
221
TYPE(Matrix_t), TARGET :: A
222
INTEGER, OPTIONAL :: ndim
223
INTEGER(KIND=AddrInt), OPTIONAL :: DotF, NormF, MatVecF, PrecF, StopcF
224
!------------------------------------------------------------------------------
225
TYPE(Matrix_t), POINTER :: Adiag,CM,PrecMat,SaveGlobalM
226
227
REAL(KIND=dp) :: dpar(HUTI_DPAR_DFLTSIZE),stopfun
228
! external stopfun
229
REAL(KIND=dp), ALLOCATABLE :: work(:,:)
230
INTEGER :: i,j,k,N,ipar(HUTI_IPAR_DFLTSIZE),wsize,istat,IterType,PCondType,ILUn,Blocks
231
LOGICAL :: Internal, NullEdges
232
LOGICAL :: ComponentwiseStopC, NormwiseStopC, RowEquilibration
233
LOGICAL :: Condition,GotIt, Refactorize,Found,GotDiagFactor,Robust
234
LOGICAL :: ComplexSystem, PseudoComplexSystem, DoFatal, LeftOriented
235
236
REAL(KIND=dp) :: ILUT_TOL, DiagFactor
237
238
TYPE(ValueList_t), POINTER :: Params
239
240
CHARACTER(:), ALLOCATABLE :: str
241
242
EXTERNAL MultigridPrec
243
EXTERNAL NormwiseBackwardError, ComponentwiseBackwardError
244
EXTERNAL NormwiseBackwardErrorGeneralized
245
EXTERNAL NormwiseBackwardError_Z
246
247
INTEGER(KIND=Addrint) :: dotProc, normProc, pcondProc, &
248
pconddProc, mvProc, iterProc, StopcProc
249
INTEGER(KIND=Addrint) :: AddrFunc
250
INTEGER :: astat
251
COMPLEX(KIND=dp), ALLOCATABLE :: xC(:), bC(:)
252
COMPLEX(KIND=dp), ALLOCATABLE :: workC(:,:)
253
EXTERNAL :: AddrFunc
254
255
INTERFACE
256
SUBROUTINE VankaCreate(A,Solver)
257
USE Types
258
TYPE(Matrix_t) :: A
259
TYPE(Solver_t) :: Solver
260
END SUBROUTINE VankaCreate
261
262
SUBROUTINE VankaPrec(u,v,ipar)
263
USE Types
264
INTEGER :: ipar(*)
265
REAL(KIND=dp) :: u(*),v(*)
266
END SUBROUTINE VankaPrec
267
268
SUBROUTINE CircuitPrec(u,v,ipar)
269
USE Types
270
INTEGER :: ipar(*)
271
REAL(KIND=dp) :: u(*),v(*)
272
END SUBROUTINE CircuitPrec
273
274
SUBROUTINE CircuitPrecComplex(u,v,ipar)
275
USE Types
276
INTEGER :: ipar(*)
277
COMPLEX(KIND=dp) :: u(*),v(*)
278
END SUBROUTINE CircuitPrecComplex
279
280
SUBROUTINE SlavePrec(u,v,ipar)
281
USE Types
282
INTEGER :: ipar(*)
283
REAL(KIND=dp) :: u(*),v(*)
284
END SUBROUTINE SlavePrec
285
286
SUBROUTINE SlavePrecComplex(u,v,ipar)
287
USE Types
288
INTEGER :: ipar(*)
289
COMPLEX(KIND=dp) :: u(*),v(*)
290
END SUBROUTINE SlavePrecComplex
291
END INTERFACE
292
!------------------------------------------------------------------------------
293
N = A % NumberOfRows
294
IF ( PRESENT(ndim) ) n=ndim
295
296
ipar = 0
297
dpar = 0.0D0
298
pconddProc = 0
299
!------------------------------------------------------------------------------
300
Params => Solver % Values
301
str = ListGetString( Params,'Linear System Iterative Method',Found )
302
IF( .NOT. Found ) THEN
303
CALL Warn('IterSolver','> Linear System Iterative Method < not found, using BiCGstab')
304
str = 'bicgstab'
305
ELSE
306
CALL Info('IterSolver','Using iterative method: '//TRIM(str),Level=9)
307
END IF
308
309
IF( ListGetLogical( Params,'Linear System Skip Complex',GotIt ) ) THEN
310
CALL Info('IterSolver','This time skipping complex treatment',Level=20)
311
A % COMPLEX = .FALSE.
312
ComplexSystem = .FALSE.
313
ELSE
314
ComplexSystem = ListGetLogical( Params,'Linear System Complex',Found )
315
IF( .NOT. Found ) ComplexSystem = A % COMPLEX
316
END IF
317
318
PseudoComplexSystem = ListGetLogical( Params,'Linear System Pseudo Complex',Found )
319
320
IF( ComplexSystem ) THEN
321
CALL Info('IterSolver','Matrix is complex valued',Level=10)
322
ELSE IF( PseudoComplexSystem ) THEN
323
CALL Info('IterSolver','Matrix is pseudo complex valued',Level=10)
324
ELSE
325
CALL Info('IterSolver','Matrix is real valued',Level=12)
326
END IF
327
328
329
SELECT CASE(str)
330
CASE('bicgstab2')
331
! NOTE:
332
! BiCGStab2 should be nearly the same as BiCGStabl with the parameter l=2, but
333
! the implementation of BiCGStabl uses the right-oriented preconditioning, while
334
! BiCGStab2 works as expected only with the left-oriented preconditioning. Due to
335
! the difference in the preconditioning the convergence histories may be quite different.
336
! The complex BiCGStab2 does not convince, use BiCGStabl instead:
337
IF (ComplexSystem ) THEN
338
IterType = ITER_BICGstabl
339
ELSE
340
IterType = ITER_BiCGStab2
341
END IF
342
CASE('bicgstabl')
343
IterType = ITER_BICGstabl
344
CASE('bicgstab')
345
IterType = ITER_BiCGStab
346
CASE('tfqmr')
347
IterType = ITER_TFQMR
348
CASE('cgs')
349
IterType = ITER_CGS
350
CASE('cg')
351
IterType = ITER_CG
352
CASE('gmres')
353
IterType = ITER_GMRES
354
CASE('sgs')
355
IterType = ITER_SGS
356
CASE('jacobi')
357
IterType = ITER_jacobi
358
CASE('richardson')
359
IterType = ITER_richardson
360
CASE('gcr')
361
IterType = ITER_GCR
362
CASE('idrs')
363
IterType = ITER_IDRS
364
CASE DEFAULT
365
IterType = ITER_BiCGStab
366
END SELECT
367
368
!------------------------------------------------------------------------------
369
370
HUTI_WRKDIM = 0
371
HUTI_PSEUDOCOMPLEX = 0
372
IF( PseudoComplexSystem ) THEN
373
HUTI_PSEUDOCOMPLEX = 1
374
IF ( ListGetLogical( Params,'Block Split Complex',Found ) ) HUTI_PSEUDOCOMPLEX = 2
375
END IF
376
Internal = .FALSE.
377
378
SELECT CASE ( IterType )
379
380
! Solvers from HUTiter
381
!-------------------------------------------------------
382
CASE (ITER_BiCGStab)
383
HUTI_WRKDIM = HUTI_BICGSTAB_WORKSIZE
384
385
CASE (ITER_BiCGStab2)
386
HUTI_WRKDIM = HUTI_BICGSTAB_2_WORKSIZE
387
388
CASE (ITER_TFQMR)
389
HUTI_WRKDIM = HUTI_TFQMR_WORKSIZE
390
391
CASE (ITER_CG)
392
HUTI_WRKDIM = HUTI_CG_WORKSIZE
393
394
CASE (ITER_CGS)
395
HUTI_WRKDIM = HUTI_CGS_WORKSIZE
396
397
CASE (ITER_GMRES)
398
HUTI_GMRES_RESTART = ListGetInteger( Params,&
399
'Linear System GMRES Restart', GotIt )
400
IF ( .NOT. GotIT ) HUTI_GMRES_RESTART = 10
401
HUTI_WRKDIM = HUTI_GMRES_WORKSIZE + HUTI_GMRES_RESTART
402
403
! Solvers from IterativeMethods.src
404
!-------------------------------------------------------
405
CASE (ITER_SGS)
406
HUTI_WRKDIM = 1
407
HUTI_SGSPARAM = ListGetConstReal( Params,'SGS Overrelaxation Factor',&
408
GotIt,minv=0.0_dp,maxv=2.0_dp)
409
IF(.NOT. GotIt) HUTI_SGSPARAM = 1.8
410
Internal = .TRUE.
411
412
CASE (ITER_Jacobi, ITER_Richardson)
413
HUTI_WRKDIM = 1
414
Internal = .TRUE.
415
416
CASE (ITER_GCR)
417
HUTI_WRKDIM = 1
418
HUTI_GCR_RESTART = ListGetInteger( Params, &
419
'Linear System GCR Restart', GotIt )
420
IF ( .NOT. GotIt ) THEN
421
i = ListGetInteger( Params,'Linear System Max Iterations', minv=1 )
422
IF( i > 200 ) THEN
423
i = 200
424
CALL Info('IterSolver','"Linear System GCR Restart" not given, setting it to '//I2S(i),Level=4)
425
END IF
426
HUTI_GCR_RESTART = i
427
END IF
428
Internal = .TRUE.
429
430
CASE (ITER_BICGSTABL)
431
HUTI_WRKDIM = 1
432
HUTI_BICGSTABL_L = ListGetInteger( Params,'BiCGstabl polynomial degree',&
433
GotIt,minv=2)
434
IF(.NOT. GotIt) HUTI_BICGSTABL_L = 2
435
Internal = .TRUE.
436
437
CASE (ITER_IDRS)
438
HUTI_WRKDIM = 1
439
HUTI_IDRS_S = ListGetInteger( Params,'IDRS parameter',GotIt,minv=1)
440
IF(.NOT. GotIt) HUTI_IDRS_S = 4
441
Internal = .TRUE.
442
443
END SELECT
444
!------------------------------------------------------------------------------
445
446
wsize = HUTI_WRKDIM
447
448
StopcProc = 0
449
IF (PRESENT(StopcF)) THEN
450
StopcProc = StopcF
451
HUTI_STOPC = HUTI_USUPPLIED_STOPC
452
ELSE
453
ComponentwiseStopC = ListGetLogical(Params,'Linear System Componentwise Backward Error',GotIt)
454
IF (ComponentwiseStopC) THEN
455
IF (ComplexSystem) THEN
456
CALL Info('IterSolver', 'Linear System Componentwise Backward Error is active')
457
CALL Fatal('IterSolver', 'Error computation does not support Linear System Complex = True')
458
END IF
459
StopcProc = AddrFunc(ComponentwiseBackwardError)
460
HUTI_STOPC = HUTI_USUPPLIED_STOPC
461
ELSE
462
NormwiseStopC = ListGetLogical(Params,'Linear System Normwise Backward Error',GotIt)
463
IF (NormwiseStopC) THEN
464
RowEquilibration = ListGetLogical(Params,'Linear System Row Equilibration',GotIt)
465
IF (RowEquilibration) THEN
466
IF (ComplexSystem) THEN
467
StopcProc = AddrFunc(NormwiseBackwardError_Z)
468
ELSE
469
StopcProc = AddrFunc(NormwiseBackwardError)
470
END IF
471
ELSE
472
IF (ComplexSystem) THEN
473
CALL Info('IterSolver', 'Linear System Normwise Backward Error is active')
474
CALL Fatal('IterSolver', 'Error computation needs Linear System Row Equilibration = True')
475
ELSE
476
StopcProc = AddrFunc(NormwiseBackwardErrorGeneralized)
477
END IF
478
END IF
479
HUTI_STOPC = HUTI_USUPPLIED_STOPC
480
ELSE
481
HUTI_STOPC = HUTI_TRESID_SCALED_BYB
482
END IF
483
END IF
484
END IF
485
HUTI_NDIM = N
486
487
HUTI_DBUGLVL = ListGetInteger( Params, &
488
'Linear System Residual Output', GotIt )
489
IF ( .NOT.Gotit ) HUTI_DBUGLVL = 1
490
491
IF ( Parenv % myPE /= 0 ) HUTI_DBUGLVL=0
492
493
HUTI_MAXIT = ListGetInteger( Params, &
494
'Linear System Max Iterations', minv=1 )
495
496
HUTI_MINIT = ListGetInteger( Params, &
497
'Linear System Min Iterations', GotIt )
498
499
IF( ComplexSystem ) THEN
500
ALLOCATE(workC(N/2,wsize), stat=istat)
501
IF ( istat /= 0 ) THEN
502
CALL Fatal( 'IterSolve', 'Memory allocation failure.' )
503
END IF
504
workC = cmplx(0,0,dp)
505
ELSE
506
ALLOCATE(work(N,wsize), stat=istat)
507
IF ( istat /= 0 ) THEN
508
CALL Fatal( 'IterSolve', 'Memory allocation failure.' )
509
END IF
510
!$OMP PARALLEL PRIVATE(j)
511
DO j=1,wsize
512
!$OMP DO
513
DO i=1,N
514
work(i,j) = real(0,dp)
515
END DO
516
!$OMP END DO
517
END DO
518
!$OMP END PARALLEL
519
END IF
520
521
IF ( (IterType == ITER_BiCGStab2 .OR. IterType == ITER_BiCGStabL .OR. &
522
IterType == ITER_BiCGStab ) .AND. ALL(x == 0.0) ) x = 1.0d-8
523
524
HUTI_INITIALX = HUTI_USERSUPPLIEDX
525
526
HUTI_TOLERANCE = ListGetCReal( Params, &
527
'Linear System Convergence Tolerance' )
528
529
HUTI_MAXTOLERANCE = ListGetCReal( Params, &
530
'Linear System Divergence Limit', GotIt)
531
IF(.NOT. GotIt) HUTI_MAXTOLERANCE = 1.0d20
532
533
IF( ListGetLogical( Params,'Linear System Robust',GotIt) ) THEN
534
HUTI_ROBUST = 1
535
HUTI_ROBUST_TOLERANCE = ListGetCReal( Params,'Linear System Robust Tolerance',GotIt)
536
IF(.NOT. GotIt ) HUTI_ROBUST_TOLERANCE = HUTI_TOLERANCE**(2.0/3.0)
537
HUTI_ROBUST_MAXTOLERANCE = ListGetCReal( Params,'Linear System Robust Limit',GotIt)
538
IF(.NOT. GotIt ) HUTI_ROBUST_MAXTOLERANCE = SQRT( HUTI_TOLERANCE )
539
HUTI_ROBUST_STEPSIZE = ListGetCReal( Params,'Linear System Robust Margin',GotIt)
540
IF(.NOT. GotIt ) HUTI_ROBUST_STEPSIZE = 1.1_dp
541
HUTI_ROBUST_MAXBADIT = ListGetInteger( Params,'Linear System Robust Max Iterations',GotIt)
542
IF(.NOT. GotIt ) HUTI_ROBUST_MAXBADIT = HUTI_MAXIT / 2
543
HUTI_ROBUST_START = ListGetInteger( Params,'Linear System Robust Start Iteration',GotIt)
544
IF(.NOT. GotIt ) HUTI_ROBUST_START = 1
545
ELSE
546
HUTI_ROBUST = 0
547
END IF
548
549
550
IF( ListGetLogical( Params,'IDRS Smoothing',GotIt) ) THEN
551
HUTI_SMOOTHING = 1
552
ELSE
553
HUTI_SMOOTHING = 0
554
END IF
555
556
557
!------------------------------------------------------------------------------
558
559
! By default the right-oriented preconditioning is applied, but BiCGStab2,
560
! GMRES and TFQMR are called with the left-oriented preconditining since
561
! the right-oriented preconditioning does not work as expected. The methods
562
! from the module IterativeMethods use always the right-oriented preconditioning:
563
!
564
SELECT CASE ( IterType )
565
CASE (ITER_BiCGStab2, ITER_GMRES, ITER_TFQMR)
566
LeftOriented = .TRUE.
567
CASE DEFAULT
568
LeftOriented = ListGetLogical(Params, 'Linear System Left Preconditioning', GotIt)
569
IF (Internal) LeftOriented = .FALSE.
570
END SELECT
571
572
573
IF ( .NOT. PRESENT(PrecF) ) THEN
574
str = ListGetString( Params, 'Linear System Preconditioning',gotit )
575
IF ( .NOT.gotit ) str = 'none'
576
577
A % Cholesky = ListGetLogical( Params,'Linear System Symmetric ILU', Gotit )
578
579
ILUn = -1
580
IF ( str == 'none' ) THEN
581
PCondType = PRECOND_NONE
582
583
ELSE IF ( str == 'diagonal' ) THEN
584
PCondType = PRECOND_DIAGONAL
585
586
ELSE IF ( str == 'ilut' ) THEN
587
ILUT_TOL = ListGetCReal( Params, &
588
'Linear System ILUT Tolerance',GotIt )
589
PCondType = PRECOND_ILUT
590
591
ELSE IF ( SEQL(str, 'ilu') ) THEN
592
ILUn = ListGetInteger( Params, &
593
'Linear System ILU Order', gotit )
594
IF ( .NOT.gotit ) THEN
595
IF(LEN(str)>=4) ILUn = ICHAR(str(4:4)) - ICHAR('0')
596
END IF
597
IF ( ILUn < 0 .OR. ILUn > 9 ) ILUn = 0
598
PCondType = PRECOND_ILUn
599
600
ELSE IF ( SEQL(str, 'bilu') ) THEN
601
ILUn = 0
602
IF(LEN(str)>=5) ILUn = ICHAR(str(5:5)) - ICHAR('0')
603
IF ( ILUn < 0 .OR. ILUn > 9 ) ILUn = 0
604
IF( Solver % Variable % Dofs == 1) THEN
605
CALL Warn('IterSolver','BILU for one dofs is equal to ILU!')
606
PCondType = PRECOND_ILUn
607
ELSE
608
PCondType = PRECOND_BILUn
609
END IF
610
611
ELSE IF ( str == 'multigrid' ) THEN
612
PCondType = PRECOND_MG
613
614
ELSE IF ( SEQL(str,'vanka') ) THEN
615
PCondType = PRECOND_VANKA
616
617
ELSE IF ( str == 'auxiliary space solver' .OR. str == 'slave' ) THEN
618
PCondType = PRECOND_SLAVE
619
620
ELSE IF ( str == 'circuit' ) THEN
621
ILUn = ListGetInteger( Params, 'Linear System ILU Order', gotit )
622
IF(.NOT.Gotit ) ILUn=-1
623
PCondType = PRECOND_Circuit
624
625
ELSE
626
PCondType = PRECOND_NONE
627
CALL Warn( 'IterSolve', 'Unknown preconditioner type: '//TRIM(str)//', feature disabled.' )
628
END IF
629
630
IF ( .NOT. ListGetLogical( Params, 'No Precondition Recompute',GotIt ) ) THEN
631
CALL ResetTimer("Prec-"//TRIM(str))
632
633
n = ListGetInteger( Params, 'Linear System Precondition Recompute', GotIt )
634
IF ( n <= 0 ) n = 1
635
636
Refactorize = ListGetLogical( Params, 'Linear System Refactorize', Gotit )
637
IF ( .NOT. Gotit ) Refactorize = .TRUE.
638
639
IF (.NOT.(ASSOCIATED(A % ILUValues).OR.ASSOCIATED(A % CILUValues)).OR. &
640
(Refactorize.AND.MOD(A % SolveCount, n)==0) ) THEN
641
642
643
IF ( A % FORMAT == MATRIX_CRS ) THEN
644
645
! Optionally one may emphasize the diagonal entries in the linear system
646
! to make the preconditioning more stable.
647
!-------------------------------------------------------------------------
648
DiagFactor = ListGetCReal( Params,'Linear System ILU Factor',GotIt )
649
GotDiagFactor = ( DiagFactor > EPSILON( DiagFactor ) )
650
IF( GotDiagFactor ) THEN
651
CALL Info('IterSolver','Applying diagonal relaxation for ILU', Level=8)
652
DiagFactor = DiagFactor + 1.0_dp
653
A % Values( A % Diag ) = DiagFactor * A % Values( A % Diag )
654
END IF
655
656
IF ( ComplexSystem ) THEN
657
IF ( PCondType == PRECOND_ILUn .OR. (PCondType==PRECOND_Circuit.AND.ILUn>=0) ) THEN
658
NullEdges = ListGetLogical(Params, 'Edge Basis', GotIt)
659
CM => A % ConstraintMatrix
660
IF(NullEdges.OR.ASSOCIATED(CM)) THEN
661
CALL Info('IterSolver','Omitting edge dofs from being target of ILUn',Level=20)
662
663
IF(ASSOCIATED(A % ILURows)) DEALLOCATE(A % ILURows)
664
IF(ASSOCIATED(A % ILUCols)) DEALLOCATE(A % ILUCols)
665
IF(ASSOCIATED(A % ILUDiag)) DEALLOCATE(A % ILUDiag)
666
IF(ASSOCIATED(A % CILUValues)) DEALLOCATE(A % CILUValues)
667
668
PrecMat => AllocateMatrix()
669
PrecMat % FORMAT = MATRIX_LIST
670
PrecMat % CIluValues => NULL()
671
672
IF(ASSOCIATED(CM)) THEN
673
DO i=CM % NumberOfRows,1,-1
674
k = i + A % NumberOfRows
675
CALL List_AddMatrixIndex( PrecMat % ListMatrix,k,k)
676
IF(MOD(k,2)==0) THEN
677
CALL List_AddMatrixIndex(PrecMat % ListMatrix, k, k-1)
678
ELSE
679
CALL List_AddMatrixIndex(PrecMat % ListMatrix, k, k+1)
680
END IF
681
682
DO j=CM % Rows(i+1)-1,CM % Rows(i),-1
683
CALL List_AddToMatrixElement( PrecMat % ListMatrix, &
684
i + A % NumberOfRows, CM % Cols(j), CM % Values(j))
685
686
CALL List_AddToMatrixElement( PrecMat % ListMatrix, &
687
CM % Cols(j), i + A % NumberOfRows, CM % Values(j))
688
END DO
689
END DO
690
END IF
691
692
k = A % NumberOfRows - A % ExtraDOFs
693
DO i=A % NumberOfRows,1,-1
694
IF(i>k) THEN
695
CALL List_AddMatrixIndex(PrecMat % ListMatrix, i, i)
696
IF(MOD(i,2)==0) THEN
697
CALL List_AddMatrixIndex(PrecMat % ListMatrix, i, i-1)
698
ELSE
699
CALL List_AddMatrixIndex(PrecMat % ListMatrix, i, i+1)
700
END IF
701
ELSE IF (NullEdges) THEN
702
CALL List_AddToMatrixElement(PrecMat % ListMatrix, i, i, 1._dp)
703
IF(MOD(i,2)==0) THEN
704
CALL List_AddMatrixIndex(PrecMat % ListMatrix, i, i-1)
705
ELSE
706
CALL List_AddMatrixIndex(PrecMat % ListMatrix, i, i+1)
707
END IF
708
END IF
709
710
DO j=A % Rows(i+1)-1,A % Rows(i),-1
711
IF (i>k .OR. A % Cols(j)>k .OR. .NOT.NullEdges) THEN
712
CALL List_AddToMatrixElement(PrecMat % ListMatrix, i, A % Cols(j), A % Values(j))
713
END IF
714
END DO
715
END DO
716
717
CALL List_ToCRSMatrix(PrecMat)
718
Condition = CRS_ComplexIncompleteLU(PrecMat,ILUn)
719
720
A % ILURows => PrecMat % IluRows
721
A % ILUCols => PrecMat % IluCols
722
A % ILUDiag => PrecMat % IluDiag
723
A % CILUvalues => PrecMat % CIluValues
724
725
DEALLOCATE(PrecMat % Values)
726
IF(.NOT.ASSOCIATED(A % ILURows,PrecMat % Rows)) DEALLOCATE(PrecMat % Rows)
727
IF(.NOT.ASSOCIATED(A % ILUCols,PrecMat % Cols)) DEALLOCATE(PrecMat % Cols)
728
IF(.NOT.ASSOCIATED(A % ILUDiag,PrecMat % Diag)) DEALLOCATE(PrecMat % Diag)
729
DEALLOCATE(PrecMat)
730
ELSE
731
Condition = CRS_ComplexIncompleteLU(A,ILUn)
732
END IF
733
ELSE IF ( PCondType == PRECOND_ILUT ) THEN
734
Condition = CRS_ComplexILUT( A,ILUT_TOL )
735
END IF
736
ELSE IF (ILUn>=0 .OR. PCondType == PRECOND_ILUT) THEN ! Not ComplexSystem
737
SELECT CASE(PCondType)
738
CASE(PRECOND_ILUn, PRECOND_Circuit)
739
NullEdges = ListGetLogical(Params, 'Edge Basis', GotIt)
740
CM => A % ConstraintMatrix
741
IF(NullEdges.OR.ASSOCIATED(CM)) THEN
742
743
IF(ASSOCIATED(A % ILURows)) DEALLOCATE(A % ILURows)
744
IF(ASSOCIATED(A % ILUCols)) DEALLOCATE(A % ILUCols)
745
IF(ASSOCIATED(A % ILUDiag)) DEALLOCATE(A % ILUDiag)
746
IF(ASSOCIATED(A % ILUValues)) DEALLOCATE(A % ILUValues)
747
748
PrecMat => AllocateMatrix()
749
PrecMat % FORMAT = MATRIX_LIST
750
751
IF(ASSOCIATED(CM)) THEN
752
DO i=CM % NumberOfRows,1,-1
753
CALL List_AddMatrixIndex( PrecMat % ListMatrix, &
754
i + A % NumberOfRows, i + A % NumberOFrows)
755
756
DO j=CM % Rows(i+1)-1,CM % Rows(i),-1
757
CALL List_AddToMatrixElement( PrecMat % ListMatrix, &
758
i + A % NumberOfRows, CM % Cols(j), CM % Values(j))
759
760
CALL List_AddToMatrixElement( PrecMat % ListMatrix, &
761
CM % Cols(j), i + A % NumberOfRows, CM % Values(j))
762
END DO
763
END DO
764
END IF
765
766
k = A % NumberOfRows - A % ExtraDOFs
767
DO i=A % NumberOfRows,1,-1
768
IF(i>k) THEN
769
CALL List_AddMatrixIndex(PrecMat % ListMatrix, i, i)
770
ELSE IF (NullEdges) THEN
771
CALL List_AddToMatrixElement(PrecMat % ListMatrix, i, i, 1._dp)
772
END IF
773
774
DO j=A % Rows(i+1)-1,A % Rows(i),-1
775
IF (i>k .OR. A % Cols(j)>k .OR. .NOT.NullEdges) THEN
776
CALL List_AddToMatrixElement(PrecMat % ListMatrix, i, A % Cols(j), A % Values(j))
777
END IF
778
END DO
779
END DO
780
781
CALL List_ToCRSMatrix(PrecMat)
782
Condition = CRS_IncompleteLU(PrecMat,ILUn,Params)
783
784
A % ILURows => PrecMat % IluRows
785
A % ILUCols => PrecMat % IluCols
786
A % ILUDiag => PrecMat % IluDiag
787
A % ILUvalues => PrecMat % IluValues
788
789
DEALLOCATE(PrecMat % Values)
790
IF(.NOT.ASSOCIATED(A % ILURows,PrecMat % Rows)) DEALLOCATE(PrecMat % Rows)
791
IF(.NOT.ASSOCIATED(A % ILUCols,PrecMat % Cols)) DEALLOCATE(PrecMat % Cols)
792
IF(.NOT.ASSOCIATED(A % ILUDiag,PrecMat % Diag)) DEALLOCATE(PrecMat % Diag)
793
DEALLOCATE(PrecMat)
794
ELSE
795
Condition = CRS_IncompleteLU(A,ILUn,Params)
796
END IF
797
CASE(PRECOND_ILUT)
798
Condition = CRS_ILUT( A,ILUT_TOL )
799
CASE(PRECOND_BILUn)
800
Blocks = Solver % Variable % Dofs
801
IF ( Blocks <= 1 ) THEN
802
Condition = CRS_IncompleteLU(A,ILUn,Params)
803
ELSE
804
IF( .NOT. ASSOCIATED( A % ILUValues ) ) THEN
805
Adiag => AllocateMatrix()
806
CALL CRS_BlockDiagonal(A,Adiag,Blocks)
807
Condition = CRS_IncompleteLU(Adiag,ILUn,Params)
808
A % ILURows => Adiag % ILURows
809
A % ILUCols => Adiag % ILUCols
810
A % ILUValues => Adiag % ILUValues
811
A % ILUDiag => Adiag % ILUDiag
812
IF (ILUn > 0) THEN
813
DEALLOCATE(Adiag % Rows,Adiag % Cols, Adiag % Diag, Adiag % Values)
814
END IF
815
DEALLOCATE( Adiag )
816
ELSE
817
Condition = CRS_IncompleteLU(A,ILUn,Params)
818
END IF
819
END IF
820
CASE(PRECOND_VANKA)
821
! CALL VankaCreate( A, SolverParam )
822
END SELECT
823
END IF
824
825
IF( GotDiagFactor ) THEN
826
CALL Info('IterSolver','Reverting diagonal relaxation for ILU', Level=10)
827
A % Values( A % Diag ) = A % Values( A % Diag ) / DiagFactor
828
END IF
829
830
ELSE
831
IF ( PCondType == PRECOND_ILUn ) THEN
832
CALL Warn( 'IterSolve', 'No ILU Preconditioner for Band Matrix format,' )
833
CALL Warn( 'IterSolve', 'using Diagonal preconditioner instead...' )
834
PCondType = PRECOND_DIAGONAL
835
END IF
836
END IF
837
END IF
838
CALL CheckTimer("Prec-"//TRIM(str),Level=8,Delete=.TRUE.)
839
END IF
840
END IF
841
842
A % SolveCount = A % SolveCount + 1
843
!------------------------------------------------------------------------------
844
845
IF ( PRESENT(MatvecF) ) THEN
846
mvProc = MatvecF
847
ELSE
848
IF ( .NOT. ComplexSystem ) THEN
849
mvProc = AddrFunc( CRS_MatrixVectorProd )
850
ELSE
851
mvProc = AddrFunc( CRS_ComplexMatrixVectorProd )
852
END IF
853
END IF
854
855
IF ( PRESENT(dotF) ) THEN
856
dotProc = dotF
857
ELSE
858
dotProc = 0
859
END IF
860
861
IF ( PRESENT(normF) ) THEN
862
normProc = normF
863
ELSE
864
normProc = 0
865
END IF
866
867
868
IF ( PRESENT(PrecF) ) THEN
869
pcondProc = PrecF
870
ELSE
871
SELECT CASE( PCondType )
872
CASE (PRECOND_NONE)
873
IF ( .NOT. ComplexSystem ) THEN
874
pcondProc = AddrFunc( pcond_dummy )
875
ELSE
876
pcondProc = AddrFunc( pcond_dummy_cmplx )
877
END IF
878
879
CASE (PRECOND_DIAGONAL)
880
IF ( .NOT. ComplexSystem ) THEN
881
pcondProc = AddrFunc( CRS_DiagPrecondition )
882
ELSE
883
pcondProc = AddrFunc( CRS_ComplexDiagPrecondition )
884
END IF
885
886
CASE (PRECOND_ILUn, PRECOND_ILUT, PRECOND_BILUn )
887
IF ( .NOT. ComplexSystem ) THEN
888
pcondProc = AddrFunc( CRS_LUPrecondition )
889
ELSE
890
pcondProc = AddrFunc( CRS_ComplexLUPrecondition )
891
END IF
892
893
CASE (PRECOND_MG)
894
pcondProc = AddrFunc( MultiGridPrec )
895
896
CASE (PRECOND_VANKA)
897
pcondProc = AddrFunc( VankaPrec )
898
899
CASE (PRECOND_Slave)
900
IF ( .NOT. ComplexSystem ) THEN
901
pcondProc = AddrFunc( SlavePrec )
902
ELSE
903
pcondProc = AddrFunc( SlavePrecComplex )
904
END IF
905
906
CASE (PRECOND_Circuit)
907
IF ( .NOT. ComplexSystem ) THEN
908
pcondProc = AddrFunc( CircuitPrec )
909
ELSE
910
pcondProc = AddrFunc( CircuitPrecComplex )
911
END IF
912
913
CASE DEFAULT
914
pcondProc = 0
915
END SELECT
916
END IF
917
918
919
IF ( .NOT. ComplexSystem ) THEN
920
SELECT CASE ( IterType )
921
922
! Solvers from HUTiter library
923
!-------------------------------------------------------
924
CASE (ITER_BiCGStab)
925
iterProc = AddrFunc( HUTI_D_BICGSTAB )
926
CASE (ITER_BiCGStab2)
927
iterProc = AddrFunc( HUTI_D_BICGSTAB_2 )
928
CASE (ITER_TFQMR)
929
iterProc = AddrFunc( HUTI_D_TFQMR )
930
CASE (ITER_CG)
931
iterProc = AddrFunc( HUTI_D_CG )
932
CASE (ITER_CGS)
933
iterProc = AddrFunc( HUTI_D_CGS )
934
CASE (ITER_GMRES)
935
iterProc = AddrFunc( HUTI_D_GMRES )
936
937
! Solvers from IterativeMethods.src
938
!-------------------------------------------------------
939
CASE (ITER_SGS)
940
iterProc = AddrFunc( itermethod_sgs )
941
CASE (ITER_JACOBI)
942
iterProc = AddrFunc( itermethod_jacobi )
943
CASE (ITER_RICHARDSON)
944
iterProc = AddrFunc( itermethod_richardson )
945
CASE (ITER_GCR)
946
iterProc = AddrFunc( itermethod_gcr )
947
CASE (ITER_BICGSTABL)
948
iterProc = AddrFunc( itermethod_bicgstabl )
949
CASE (ITER_IDRS)
950
iterProc = AddrFunc( itermethod_idrs )
951
952
END SELECT
953
954
IF( Internal ) THEN
955
956
IF( PseudoComplexSystem ) THEN
957
IF( HUTI_PSEUDOCOMPLEX == 1 ) THEN
958
CALL Info('IterSolver','Setting dot product function to: PseudoZDotProd',Level=15)
959
dotProc = AddrFunc( PseudoZDotProd )
960
ELSE
961
CALL Info('IterSolver','Setting dot product function to: PseudoZDotProd2',Level=15)
962
dotProc = AddrFunc( PseudoZDotProd2 )
963
END IF
964
ELSE
965
IF ( dotProc == 0 ) dotProc = AddrFunc(ddot)
966
END IF
967
IF ( normProc == 0 ) normproc = AddrFunc(dnrm2)
968
IF( HUTI_DBUGLVL == 0) HUTI_DBUGLVL = HUGE( HUTI_DBUGLVL )
969
END IF
970
971
ELSE
972
HUTI_NDIM = HUTI_NDIM / 2
973
SELECT CASE ( IterType )
974
975
! Solvers from HUTiter library
976
!-------------------------------------------------------
977
CASE (ITER_BiCGStab)
978
iterProc = AddrFunc( HUTI_Z_BICGSTAB )
979
CASE (ITER_BiCGStab2)
980
iterProc = AddrFunc( HUTI_Z_BICGSTAB_2 )
981
CASE (ITER_TFQMR)
982
iterProc = AddrFunc( HUTI_Z_TFQMR )
983
CASE (ITER_CG)
984
iterProc = AddrFunc( HUTI_Z_CG )
985
CASE (ITER_CGS)
986
iterProc = AddrFunc( HUTI_Z_CGS )
987
CASE (ITER_GMRES)
988
iterProc = AddrFunc( HUTI_Z_GMRES )
989
990
! Solvers from IterativeMethods.src
991
!-------------------------------------------------------
992
CASE (ITER_GCR)
993
iterProc = AddrFunc( itermethod_z_gcr )
994
CASE (ITER_BICGSTABL)
995
iterProc = AddrFunc( itermethod_z_bicgstabl )
996
CASE (ITER_IDRS)
997
iterProc = AddrFunc( itermethod_z_idrs )
998
CASE DEFAULT
999
CALL Fatal('IterSolver', 'Complex arithmetic version of the given linear solver is not available')
1000
END SELECT
1001
1002
IF( Internal ) THEN
1003
IF ( dotProc == 0 ) dotProc = AddrFunc(zdotc)
1004
IF ( normProc == 0 ) normproc = AddrFunc(dznrm2)
1005
IF( HUTI_DBUGLVL == 0) HUTI_DBUGLVL = HUGE( HUTI_DBUGLVL )
1006
END IF
1007
1008
END IF
1009
1010
!------------------------------------------------------------------------------
1011
1012
stack_pos = stack_pos+1
1013
IF(stack_pos>stack_max) THEN
1014
CALL Fatal('IterSolver', 'Recursion too deep ('//I2S(stack_pos)//' vs '//I2S(stack_max)//')')
1015
ELSE IF(stack_pos<=0) THEN
1016
CALL Fatal('IterSolver', 'eh')
1017
END IF
1018
FirstCall(stack_pos) = .TRUE.
1019
1020
SaveGlobalM => GlobalMatrix
1021
GlobalMatrix => A
1022
1023
IF ( ComplexSystem ) THEN
1024
! Associate xC and bC with complex variables
1025
ALLOCATE(xC(HUTI_NDIM), bC(HUTI_NDIM), STAT=astat)
1026
IF (astat /= 0) THEN
1027
CALL Fatal('IterSolve','Unable to allocate memory for complex arrays')
1028
END IF
1029
! Initialize xC and bC
1030
DO i=1,HUTI_NDIM
1031
xC(i) = cmplx(x(2*i-1),x(2*i),dp)
1032
END DO
1033
DO i=1,HUTI_NDIM
1034
bC(i) = cmplx(b(2*i-1),b(2*i),dp)
1035
END DO
1036
1037
CALL Info('IterSolver','Calling complex iterative solver',Level=32)
1038
1039
IF (LeftOriented) THEN
1040
CALL IterCall( iterProc, xC, bC, ipar, dpar, workC, &
1041
mvProc, pcondProc, pconddProc, dotProc, normProc, stopcProc )
1042
ELSE
1043
CALL IterCall( iterProc, xC, bC, ipar, dpar, workC, &
1044
mvProc, pconddProc, pcondProc, dotProc, normProc, stopcProc )
1045
END IF
1046
1047
! Copy result back
1048
DO i=1,HUTI_NDIM
1049
x(2*i-1) = REAL(REAL(xC(i)),dp)
1050
x(2*i) = REAL(AIMAG(xC(i)),dp)
1051
END DO
1052
DEALLOCATE(bC,xC)
1053
ELSE
1054
CALL Info('IterSolver','Calling real-valued iterative solver',Level=32)
1055
1056
IF (LeftOriented) THEN
1057
CALL IterCall( iterProc, x, b, ipar, dpar, work, &
1058
mvProc, pcondProc, pconddProc, dotProc, normProc, stopcProc )
1059
ELSE
1060
CALL IterCall( iterProc, x, b, ipar, dpar, work, &
1061
mvProc, pconddProc, pcondProc, dotProc, normProc, stopcProc )
1062
END IF
1063
ENDIF
1064
1065
GlobalMatrix => SaveGlobalM
1066
1067
stack_pos=stack_pos-1
1068
1069
IF ( ComplexSystem ) HUTI_NDIM = HUTI_NDIM * 2
1070
1071
!------------------------------------------------------------------------------
1072
IF ( HUTI_INFO == HUTI_CONVERGENCE ) THEN
1073
IF( ASSOCIATED( Solver % Variable ) ) THEN
1074
Solver % Variable % LinConverged = 1
1075
END IF
1076
ELSE
1077
CALL Info('IterSolve','Returned return code: '//I2S(HUTI_INFO),Level=15)
1078
IF( HUTI_INFO == HUTI_DIVERGENCE ) THEN
1079
CALL NumericalError( 'IterSolve', 'System diverged over maximum tolerance.')
1080
ELSE IF( HUTI_INFO == HUTI_MAXITER ) THEN
1081
DoFatal = ListGetLogical( Params,'Linear System Abort Not Converged',Found )
1082
IF(.NOT. Found ) DoFatal = .TRUE.
1083
IF( DoFatal ) THEN
1084
CALL NumericalError('IterSolve','Too many iterations were needed.')
1085
ELSE
1086
CALL Info('IterSolve','Linear iteration did not converge to tolerance',Level=6)
1087
END IF
1088
ELSE IF( HUTI_INFO == HUTI_HALTED ) THEN
1089
CALL Warn('IterSolve','Iteration halted due to problem in algorithm, trying to continue')
1090
END IF
1091
IF( ASSOCIATED( Solver % Variable ) ) THEN
1092
Solver % Variable % LinConverged = 0
1093
END IF
1094
END IF
1095
!------------------------------------------------------------------------------
1096
IF ( ComplexSystem ) THEN
1097
DEALLOCATE( workC )
1098
ELSE
1099
DEALLOCATE( work )
1100
END IF
1101
1102
!------------------------------------------------------------------------------
1103
END SUBROUTINE IterSolver
1104
!------------------------------------------------------------------------------
1105
1106
!-----------------------------------------------------------------------
1107
!> This routine may be used to either inform user or terminate following
1108
!> convergence/numerical issues, based on a flag in the SIF. Default
1109
!> behaviour terminates execution.
1110
!-----------------------------------------------------------------------
1111
SUBROUTINE NumericalError( Caller, String, IsFatal )
1112
!-----------------------------------------------------------------------
1113
CHARACTER(LEN=*) :: Caller, String
1114
LOGICAL, OPTIONAL :: IsFatal
1115
!-----------------------------------------------------------------------
1116
LOGICAL :: DoFatal, Found
1117
!-----------------------------------------------------------------------
1118
1119
!Fatality logic:
1120
! 1) Respect calling routine's wishes if present
1121
! 2) Respect solver specific option if present
1122
! 3) Respect global abort flag if present
1123
! 4) Otherwise fatal (backwards compatibility)
1124
1125
IF(PRESENT(IsFatal)) THEN
1126
DoFatal = IsFatal
1127
ELSE
1128
DoFatal = ListGetLogical(CurrentModel % Simulation,&
1129
'Global Abort Not Converged',Found)
1130
IF(.NOT. Found ) DoFatal = .TRUE.
1131
END IF
1132
1133
IF(DoFatal) THEN
1134
CALL Fatal(Caller,'Numerical Error: '//TRIM(String))
1135
ELSE
1136
CALL Warn(Caller,'Numerical Error: '//TRIM(String))
1137
END IF
1138
1139
!-----------------------------------------------------------------------
1140
END SUBROUTINE NumericalError
1141
!-----------------------------------------------------------------------
1142
1143
1144
END MODULE IterSolve
1145
1146
!> \} ElmerLib
1147
1148