/
SolverUtils.F90
22996 lines (18984 loc) · 818 KB
/
SolverUtils.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
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!/*****************************************************************************/
! *
! * Elmer, A Finite Element Software for Multiphysical Problems
! *
! * Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland
! *
! * This library is free software; you can redistribute it and/or
! * modify it under the terms of the GNU Lesser General Public
! * License as published by the Free Software Foundation; either
! * version 2.1 of the License, or (at your option) any later version.
! *
! * This library is distributed in the hope that it will be useful,
! * but WITHOUT ANY WARRANTY; without even the implied warranty of
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
! * Lesser General Public License for more details.
! *
! * You should have received a copy of the GNU Lesser General Public
! * License along with this library (in file ../LGPL-2.1); if not, write
! * to the Free Software Foundation, Inc., 51 Franklin Street,
! * Fifth Floor, Boston, MA 02110-1301 USA
! *
! *****************************************************************************/
!
!/******************************************************************************
! *
! * Utilities for *Solver - routines
! *
! ******************************************************************************
! *
! * Authors: Juha Ruokolainen
! * Email: Juha.Ruokolainen@csc.fi
! * Web: http://www.csc.fi/elmer
! * Address: CSC - IT Center for Science Ltd.
! * Keilaranta 14
! * 02101 Espoo, Finland
! *
! * Original Date: 28 Sep 1998
! *
! *****************************************************************************/
!> Basic utilities used by individual solvers.
!------------------------------------------------------------------------------
!> \ingroup ElmerLib
!> \{
MODULE SolverUtils
#include "../config.h"
USE LoadMod
USE DirectSolve
USE Multigrid
USE IterSolve
USE ElementUtils
USE ComponentUtils
USE TimeIntegrate
USE ModelDescription
USE MeshUtils
USE ParallelUtils
USE ParallelEigenSolve
USE ListMatrix
USE CRSMatrix
USE MatrixAssembly
IMPLICIT NONE
CONTAINS
!> Initialize matrix structure and vector to zero initial value.
!------------------------------------------------------------------------------
SUBROUTINE InitializeToZero( A, ForceVector )
!------------------------------------------------------------------------------
TYPE(Matrix_t), POINTER :: A !< Matrix to be initialized
REAL(KIND=dp) :: ForceVector(:) !< vector to be initialized
!------------------------------------------------------------------------------
INTEGER :: i,dim
LOGICAL :: Found, AnyNT, AnyProj, DoDisplaceMesh
TYPE(Solver_t), POINTER :: Solver
TYPE(NormalTangential_t), POINTER :: NT
CHARACTER(LEN=MAX_NAME_LEN) :: str
!------------------------------------------------------------------------------
CALL Info('InitializeToZero','Initializing the linear system to zero',Level=12)
IF ( ASSOCIATED( A ) ) THEN
SELECT CASE( A % FORMAT )
CASE( MATRIX_CRS )
CALL CRS_ZeroMatrix( A )
CASE( MATRIX_BAND,MATRIX_SBAND )
CALL Band_ZeroMatrix( A )
END SELECT
IF ( ASSOCIATED(A % PrecValues) ) THEN
A % PrecValues(:) = 0._dp
END IF
IF ( ASSOCIATED( A % MassValues ) ) THEN
A % MassValues(:) = 0.d0
END IF
IF ( ASSOCIATED( A % DampValues ) ) THEN
A % DampValues(:) = 0.d0
END IF
IF ( ASSOCIATED( A % Force ) ) THEN
A % Force(:,1) = 0.0d0
END IF
IF ( ASSOCIATED( A % RHS_im ) ) THEN
A % RHS_im(:) = 0.0d0
END IF
END IF
ForceVector = 0.0d0
Solver => CurrentModel % Solver
IF ( Solver % Variable % DOFs <= 1 ) RETURN
str = 'Normal-Tangential'
IF ( SEQL(Solver % Variable % Name, 'flow solution') ) THEN
str = TRIM(str) // ' Velocity'
ELSE
str = TRIM(str) // ' ' // GetVarName(Solver % Variable)
END IF
AnyNT = ListGetLogicalAnyBC( CurrentModel, str )
AnyProj = ListGetLogicalAnyBC( CurrentModel, 'Mortar BC Nonlinear')
IF( .NOT. (AnyNT .OR. AnyProj ) ) RETURN
DoDisplaceMesh = ListGetLogical( Solver % Values,'Displace Mesh At Init',Found )
IF( DoDisplaceMesh ) THEN
CALL Info('InitializeToZero','Displacing mesh for nonlinear projectors',Level=8)
CALL DisplaceMesh( Solver % Mesh, Solver % variable % Values, 1, &
Solver % Variable % Perm, Solver % variable % Dofs )
END IF
IF( AnyNT ) THEN
dim = CoordinateSystemDimension()
NT => CurrentModel % Solver % NormalTangential
NT % NormalTangentialNOFNodes = 0
NT % NormalTangentialName = str
CALL CheckNormalTangentialBoundary( CurrentModel, NT % NormalTangentialName, &
NT % NormalTangentialNOFNodes, NT % BoundaryReorder, &
NT % BoundaryNormals, NT % BoundaryTangent1, NT % BoundaryTangent2, dim )
CALL AverageBoundaryNormals( CurrentModel, NT % NormalTangentialName, &
NT % NormalTangentialNOFNodes, NT % BoundaryReorder, &
NT % BoundaryNormals, NT % BoundaryTangent1, NT % BoundaryTangent2, &
dim )
END IF
IF( AnyProj ) THEN
CALL GenerateProjectors(CurrentModel,Solver,Nonlinear = .TRUE. )
END IF
IF( DoDisplaceMesh ) THEN
CALL DisplaceMesh( Solver % Mesh, Solver % variable % Values, -1, &
Solver % Variable % Perm, Solver % variable % Dofs )
END IF
!------------------------------------------------------------------------------
END SUBROUTINE InitializeToZero
!------------------------------------------------------------------------------
!> Matrix vector multiplication of sparse matrices.
!------------------------------------------------------------------------------
SUBROUTINE MatrixVectorMultiply( A,u,v )
!------------------------------------------------------------------------------
TYPE(Matrix_t) :: A
INTEGER :: n
REAL(KIND=dp), DIMENSION(:) CONTIG :: u,v
!------------------------------------------------------------------------------
SELECT CASE( A % FORMAT )
CASE( MATRIX_CRS )
CALL CRS_MatrixVectorMultiply( A,u,v )
CASE( MATRIX_BAND,MATRIX_SBAND )
CALL Band_MatrixVectorMultiply( A,u,v )
CASE( MATRIX_LIST )
CALL Warn('MatrixVectorMultiply','Not implemented for List matrix type')
END SELECT
!------------------------------------------------------------------------------
END SUBROUTINE MatrixVectorMultiply
!------------------------------------------------------------------------------
!> Matrix vector multiplication of sparse matrices.
!------------------------------------------------------------------------------
SUBROUTINE MaskedMatrixVectorMultiply( A,u,v,ActiveRow,ActiveCol )
!------------------------------------------------------------------------------
TYPE(Matrix_t) :: A
INTEGER :: n
REAL(KIND=dp), DIMENSION(:) CONTIG :: u,v
LOGICAL, DIMENSION(:) :: ActiveRow
LOGICAL, DIMENSION(:) :: ActiveCol
!------------------------------------------------------------------------------
SELECT CASE( A % FORMAT )
CASE( MATRIX_CRS )
CALL CRS_MaskedMatrixVectorMultiply( A,u,v,ActiveRow, ActiveCol )
CASE DEFAULT
CALL Fatal('MaskedMatrixVectorMultiply','Not implemented for List matrix type')
END SELECT
!------------------------------------------------------------------------------
END SUBROUTINE MaskedMatrixVectorMultiply
!------------------------------------------------------------------------------
!> Matrix vector multiplication of sparse matrices.
!------------------------------------------------------------------------------
SUBROUTINE TransposeMatrixVectorMultiply( A,u,v )
!------------------------------------------------------------------------------
TYPE(Matrix_t) :: A
INTEGER :: n
REAL(KIND=dp), DIMENSION(:) CONTIG :: u,v
!------------------------------------------------------------------------------
SELECT CASE( A % FORMAT )
CASE( MATRIX_CRS )
CALL CRS_TransposeMatrixVectorMultiply( A,u,v )
CASE DEFAULT
CALL Fatal('TransposeMatrixVectorMultiply','Not implemented for other than CRS type')
END SELECT
!------------------------------------------------------------------------------
END SUBROUTINE TransposeMatrixVectorMultiply
!------------------------------------------------------------------------------
!> Search faces between passive / non-passive domains; add to boundary
!> elements with given bc-id.
!------------------------------------------------------------------------------
SUBROUTINE GetPassiveBoundary(Model,Mesh,BcId)
!------------------------------------------------------------------------------
TYPE(Model_t) :: Model
INTEGER :: BcId
TYPE(Mesh_t) :: Mesh
INTEGER, ALLOCATABLE :: arr(:)
INTEGER :: i,j,n,cnt,ind, sz
LOGICAL :: L1,L2
TYPE(Element_t), POINTER :: Faces(:), Telems(:), Face, P1, P2
CALL FindMeshEdges(Mesh,.FALSE.)
SELECT CASE(Mesh % MeshDim)
CASE(2)
Faces => Mesh % Edges
n = Mesh % NumberOfEdges
CASE(3)
Faces => Mesh % Faces
n = Mesh % NumberOfFaces
END SELECT
ALLOCATE(arr(n)); cnt=0
DO i=1,n
P1 => Faces(i) % BoundaryInfo % Right
P2 => Faces(i) % BoundaryInfo % Left
IF ( .NOT. ASSOCIATED(P1) .OR. .NOT. ASSOCIATED(P2) ) CYCLE
L1 = CheckPassiveElement(P1)
L2 = CheckPassiveElement(P2)
IF ( L1.NEQV.L2) THEN
cnt = cnt+1
arr(cnt) = i
END IF
END DO
sz = Mesh % NumberOfBulkElements + Mesh % NumberOFBoundaryElements - &
Mesh % PassBCcnt
IF ( sz+cnt>SIZE(Mesh % Elements) ) THEN
Telems => Mesh % Elements
ALLOCATE(Mesh % Elements(sz+cnt))
IF ( ASSOCIATED(Model % Elements,Telems) ) &
Model % Elements => Mesh % Elements
Mesh % Elements(1:sz) = Telems
! fix boundary element parent pointers to use new array ...
! --------------------------------------------------------
DO i=1,Mesh % NumberOfBoundaryElements-Mesh % PassBCcnt
ind = i+Mesh % NumberOfBulkElements
Face => Mesh % Elements(ind)
IF ( ASSOCIATED(Face % BoundaryInfo % Left) ) &
Face % BoundaryInfo % Left => &
Mesh % Elements(Face % BoundaryInfo % Left % ElementIndex)
IF ( ASSOCIATED(Face % BoundaryInfo % Right ) ) &
Face % BoundaryInfo % Right => &
Mesh % Elements(Face % BoundaryInfo % Right % ElementIndex)
END DO
! ...likewise for faces (edges).
! -------------------------------
DO i=1,n
Face => Faces(i)
IF ( ASSOCIATED(Face % BoundaryInfo % Left) ) &
Face % BoundaryInfo % Left => &
Mesh % Elements(Face % BoundaryInfo % Left % ElementIndex)
IF ( ASSOCIATED(Face % BoundaryInfo % Right ) ) &
Face % BoundaryInfo % Right => &
Mesh % Elements(Face % BoundaryInfo % Right % ElementIndex)
END DO
DEALLOCATE(Telems)
END IF
DO i=1,cnt
sz = sz+1
Mesh % Elements(sz) = Faces(arr(i))
Mesh % Elements(sz) % Copy = .TRUE.
Mesh % Elements(sz) % ElementIndex = sz
Mesh % Elements(sz) % BoundaryInfo % Constraint = BcId
END DO
Mesh % NumberOfBoundaryElements = Mesh % NumberOfBoundaryElements - &
Mesh % PassBCcnt + cnt
Mesh % PassBCcnt = cnt
IF ( ASSOCIATED(Model % Elements,Mesh % Elements) ) &
Model % NumberOfBoundaryElements = Mesh % NumberOfBoundaryElements
!------------------------------------------------------------------------------
END SUBROUTINE GetPassiveBoundary
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
!> For time dependent simulations add the time derivative coefficient terms
!> to the local matrix containing other coefficients.
!------------------------------------------------------------------------------
SUBROUTINE Add1stOrderTime( MassMatrix, StiffMatrix, &
Force, dt, n, DOFs, NodeIndexes, Solver, UElement )
!------------------------------------------------------------------------------
REAL(KIND=dp) :: MassMatrix(:,:) !< Local mass matrix.
REAL(KIND=dp) :: StiffMatrix(:,:) !< Local stiffness matrix.
REAL(KIND=dp) :: Force(:) !< Local right-hand-side vector.
REAL(KIND=dp) :: dt !< Simulation timestep size
INTEGER :: n !< number of element nodes
INTEGER :: DOFs !< variable degrees of freedom
INTEGER :: NodeIndexes(:) !< element nodes
TYPE(Solver_t) :: Solver !< Solver structure.
TYPE(Element_t), TARGET, OPTIONAL :: UElement !< Element structure
!------------------------------------------------------------------------------
LOGICAL :: GotIt
INTEGER :: i,j,k,l,m,Order
REAL(KIND=dp) :: s, t, zeta
CHARACTER(LEN=MAX_NAME_LEN) :: Method
REAL(KIND=dp) :: PrevSol(DOFs*n,Solver % Order), CurSol(DOFs*n), LForce(n*DOFs)
TYPE(Variable_t), POINTER :: DtVar
REAL(KIND=dp) :: Dts(Solver % Order)
LOGICAL :: ConstantDt
TYPE(Element_t), POINTER :: Element
!------------------------------------------------------------------------------
INTEGER :: PredCorrOrder !< Order of predictor-corrector scheme
IF ( PRESENT(UElement) ) THEN
Element => UElement
ELSE
Element => CurrentModel % CurrentElement
END IF
IF ( Solver % Matrix % Lumped ) THEN
#ifndef OLD_LUMPING
s = 0.d0
t = 0.d0
DO i=1,n*DOFs
DO j=1,n*DOFs
s = s + MassMatrix(i,j)
IF (i /= j) THEN
MassMatrix(i,j) = 0.d0
END IF
END DO
t = t + MassMatrix(i,i)
END DO
DO i=1,n
DO j=1,DOFs
K = DOFs * (i-1) + j
L = DOFs * (NodeIndexes(i)-1) + j
IF ( t /= 0.d0 ) THEN
MassMatrix(K,K) = MassMatrix(K,K) * s / t
END IF
END DO
END DO
#else
DO i=1,n*DOFs
s = 0.0d0
DO j = 1,n*DOFs
s = s + MassMatrix(i,j)
MassMatrix(i,j) = 0.0d0
END DO
MassMatrix(i,i) = s
END DO
DO i=1,n
DO j=1,DOFs
K = DOFs * (i-1) + j
L = DOFs * (NodeIndexes(i)-1) + j
END DO
END DO
#endif
END IF
!------------------------------------------------------------------------------
Order = MIN(Solver % DoneTime, Solver % Order)
DO i=1,n
DO j=1,DOFs
K = DOFs * (i-1) + j
L = DOFs * (NodeIndexes(i)-1) + j
DO m=1, Order
PrevSol(K,m) = Solver % Variable % PrevValues(L,m)
END DO
CurSol(K) = Solver % Variable % Values(L)
END DO
END DO
LForce(1:n*DOFs) = Force(1:n*DOFs)
CALL UpdateGlobalForce( Solver % Matrix % Force(:,1), LForce, &
n, DOFs, NodeIndexes, UElement=Element )
!------------------------------------------------------------------------------
!PrevSol(:,Order) needed for BDF
Method = ListGetString( Solver % Values, 'Timestepping Method', GotIt )
SELECT CASE( Method )
CASE( 'fs' )
CALL FractionalStep( n*DOFs, dt, MassMatrix, StiffMatrix, Force, &
PrevSol(:,1), Solver % Beta, Solver )
CASE('bdf')
Dts(1) = Dt
ConstantDt = .TRUE.
IF(Order > 1) THEN
DtVar => VariableGet( Solver % Mesh % Variables, 'Timestep size' )
DO i=2,Order
Dts(i) = DtVar % PrevValues(1,i-1)
IF(ABS(Dts(i)-Dts(1)) > 1.0d-6 * Dts(1)) ConstantDt = .FALSE.
END DO
END IF
IF(ConstantDt) THEN
CALL BDFLocal( n*DOFs, dt, MassMatrix, StiffMatrix, Force, PrevSol, &
Order )
ELSE
CALL VBDFLocal( n*DOFs, dts, MassMatrix, StiffMatrix, Force, PrevSol, &
Order )
END IF
CASE('runge-kutta')
CALL RungeKutta( n*DOFs, dt, MassMatrix, StiffMatrix, Force, &
PrevSol(:,1), CurSol )
CASE('adams-bashforth')
zeta = ListGetConstReal( Solver % Values, 'Adams Zeta', GotIt )
IF ( .NOT. Gotit) zeta = 1.0_dp
PredCorrOrder = ListGetInteger( Solver % Values, &
'Predictor-Corrector Scheme Order', GotIt)
IF (.NOT. GotIt) PredCorrOrder = 2
PredCorrOrder = MIN(PredCorrOrder, Solver % DoneTime /2)
CALL AdamsBashforth( n*DOFs, dt, MassMatrix, StiffMatrix, Force, &
PrevSol(:,1), zeta, PredCorrOrder)
CASE('adams-moulton')
PredCorrOrder = ListGetInteger( Solver % Values, &
'Predictor-Corrector Scheme Order', GotIt)
IF (.NOT. GotIt) PredCorrOrder = 2
PredCorrOrder = MIN(PredCorrOrder, Solver % DoneTime /2)
CALL AdamsMoulton( n*DOFs, dt, MassMatrix, StiffMatrix, Force, &
PrevSol, PredCorrOrder )
CASE DEFAULT
CALL NewmarkBeta( n*DOFs, dt, MassMatrix, StiffMatrix, Force, &
PrevSol(:,1), Solver % Beta )
END SELECT
!------------------------------------------------------------------------------
END SUBROUTINE Add1stOrderTime
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
!> For time dependent simulations add the time derivative coefficient terms
!> to the global matrix containing other coefficients.
!------------------------------------------------------------------------------
SUBROUTINE Add1stOrderTime_CRS( Matrix, Force, dt, Solver )
!------------------------------------------------------------------------------
TYPE(Matrix_t), POINTER :: Matrix !< Global matrix (including stiffness and mass)
REAL(KIND=dp) :: Force(:) !< Global right-hand-side vector.
REAL(KIND=dp) :: dt !< Simulation timestep size
TYPE(Solver_t) :: Solver !< Solver structure.
!------------------------------------------------------------------------------
LOGICAL :: GotIt
INTEGER :: i,j,k,l,m,n,Order
REAL(KIND=dp) :: s, t, msum
CHARACTER(LEN=MAX_NAME_LEN) :: Method
TYPE(Variable_t), POINTER :: DtVar
REAL(KIND=dp) :: Dts(Solver % Order)
REAL(KIND=dp), POINTER :: PrevSol(:,:), ML(:), CurrSol(:)
INTEGER, POINTER :: Rows(:), Cols(:)
LOGICAL :: ConstantDt, Lumped, Found
!------------------------------------------------------------------------------
CALL Info('Add1stOrderTime_CRS','Adding time discretization to CRS matrix',Level=20)
!------------------------------------------------------------------------------
Order = MIN(Solver % DoneTime, Solver % Order)
Method = ListGetString( Solver % Values, 'Timestepping Method', GotIt )
CurrSol => Solver % Variable % Values
PrevSol => Solver % Variable % PrevValues
SELECT CASE( Method )
CASE( 'fs' )
CALL FractionalStep_CRS( dt, Matrix, Force, PrevSol(:,1), Solver )
CASE('bdf')
ConstantDt = .TRUE.
IF(Order > 1) THEN
Dts(1) = Dt
DtVar => VariableGet( Solver % Mesh % Variables, 'Timestep size' )
DO i=2,Order
Dts(i) = DtVar % PrevValues(1,i-1)
IF(ABS(Dts(i)-Dts(1)) > 1.0d-6 * Dts(1)) ConstantDt = .FALSE.
END DO
END IF
IF(ConstantDt) THEN
CALL BDF_CRS( dt, Matrix, Force, PrevSol, Order )
ELSE
CALL VBDF_CRS( dts, Matrix, Force, PrevSol, Order )
END IF
CASE('runge-kutta')
CALL RungeKutta_CRS( dt, Matrix, Force, PrevSol(:,1), CurrSol )
CASE DEFAULT
CALL NewmarkBeta_CRS( dt, Matrix, Force, PrevSol(:,1), &
Solver % Beta )
END SELECT
!------------------------------------------------------------------------------
END SUBROUTINE Add1stOrderTime_CRS
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
!> For time dependent simulations add the time derivative coefficient terms
!> to the matrix containing other coefficients.
!------------------------------------------------------------------------------
SUBROUTINE Add2ndOrderTime( MassMatrix, DampMatrix, StiffMatrix, &
Force, dt, n, DOFs, NodeIndexes, Solver )
!------------------------------------------------------------------------------
REAL(KIND=dp) :: MassMatrix(:,:) !< Local mass matrix.
REAL(KIND=dp) :: DampMatrix(:,:) !< Local damping matrix.
REAL(KIND=dp) :: StiffMatrix(:,:) !< Local stiffness matrix.
REAL(KIND=dp) :: Force(:) !< Local right-hand-side vector.
REAL(KIND=dp) :: dt !< Simulation timestep size
INTEGER :: n !< number of element nodes
INTEGER :: DOFs !< variable degrees of freedom
INTEGER :: NodeIndexes(:) !< element nodes
TYPE(Solver_t) :: Solver !< Solver structure.
!------------------------------------------------------------------------------
LOGICAL :: GotIt
INTEGER :: i,j,k,l
CHARACTER(LEN=MAX_NAME_LEN) :: Method
REAL(KIND=dp) :: s,t
REAL(KIND=dp) :: X(DOFs*n),V(DOFs*N),A(DOFs*N),LForce(n*DOFs)
!------------------------------------------------------------------------------
IF ( Solver % Matrix % Lumped ) THEN
!------------------------------------------------------------------------------
#ifndef OLD_LUMPING
s = 0.d0
t = 0.d0
DO i=1,n*DOFs
DO j=1,n*DOFs
s = s + MassMatrix(i,j)
IF (i /= j) THEN
MassMatrix(i,j) = 0.d0
END IF
END DO
t = t + MassMatrix(i,i)
END DO
DO i=1,n
DO j=1,DOFs
K = DOFs * (i-1) + j
IF ( t /= 0.d0 ) THEN
MassMatrix(K,K) = MassMatrix(K,K) * s / t
END IF
END DO
END DO
s = 0.d0
t = 0.d0
DO i=1,n*DOFs
DO j=1,n*DOFs
s = s + DampMatrix(i,j)
IF (i /= j) THEN
DampMatrix(i,j) = 0.d0
END IF
END DO
t = t + DampMatrix(i,i)
END DO
DO i=1,n
DO j=1,DOFs
K = DOFs * (i-1) + j
IF ( t /= 0.d0 ) THEN
DampMatrix(K,K) = DampMatrix(K,K) * s / t
END IF
END DO
END DO
#else
!------------------------------------------------------------------------------
! Lump the second order time derivative terms ...
!------------------------------------------------------------------------------
DO i=1,n*DOFs
s = 0.0D0
DO j=1,n*DOFs
s = s + MassMatrix(i,j)
MassMatrix(i,j) = 0.0d0
END DO
MassMatrix(i,i) = s
END DO
!------------------------------------------------------------------------------
! ... and the first order terms.
!------------------------------------------------------------------------------
DO i=1,n*DOFs
s = 0.0D0
DO j=1,n*DOFs
s = s + DampMatrix(i,j)
DampMatrix(i,j) = 0.0d0
END DO
DampMatrix(i,i) = s
END DO
#endif
!------------------------------------------------------------------------------
END IF
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
! Get previous solution vectors and update current force
!-----------------------------------------------------------------------------
DO i=1,n
DO j=1,DOFs
K = DOFs * (i-1) + j
IF ( NodeIndexes(i) > 0 ) THEN
L = DOFs * (NodeIndexes(i)-1) + j
SELECT CASE(Method)
CASE DEFAULT
X(K) = Solver % Variable % PrevValues(L,3)
V(K) = Solver % Variable % PrevValues(L,4)
A(K) = Solver % Variable % PrevValues(L,5)
END SELECT
END IF
END DO
END DO
LForce(1:n*DOFs) = Force(1:n*DOFs)
CALL UpdateGlobalForce( Solver % Matrix % Force(:,1), LForce, &
n, DOFs, NodeIndexes )
!------------------------------------------------------------------------------
Method = ListGetString( Solver % Values, 'Timestepping Method', GotIt )
SELECT CASE(Method)
CASE DEFAULT
CALL Bossak2ndOrder( n*DOFs, dt, MassMatrix, DampMatrix, StiffMatrix, &
Force, X, V, A, Solver % Alpha )
END SELECT
!------------------------------------------------------------------------------
END SUBROUTINE Add2ndOrderTime
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
!> Update the right-hand-side of the global equation by adding the local entry.
!------------------------------------------------------------------------------
SUBROUTINE UpdateTimeForce( StiffMatrix, &
ForceVector, LocalForce, n, NDOFs, NodeIndexes )
!------------------------------------------------------------------------------
TYPE(Matrix_t), POINTER :: StiffMatrix !< Global stiffness matrix.
REAL(KIND=dp) :: LocalForce(:) !< Local right-hand-side vector.
REAL(KIND=dp) :: ForceVector(:) !< Global right-hand-side vector.
INTEGER :: n !< number of element nodes
INTEGER :: nDOFs !< variable degrees of freedom
INTEGER :: NodeIndexes(:) !< Element node to global node numbering mapping.
!------------------------------------------------------------------------------
INTEGER :: i,j,k
!------------------------------------------------------------------------------
CALL UpdateGlobalForce( StiffMatrix % Force(:,1), LocalForce, &
n, NDOFs, NodeIndexes )
LocalForce = 0.0d0
!------------------------------------------------------------------------------
END SUBROUTINE UpdateTimeForce
!------------------------------------------------------------------------------
!------------------------------------------------------------------------------
!> Add element local matrices & vectors to global matrices and vectors.
!------------------------------------------------------------------------------
SUBROUTINE UpdateGlobalEquations( StiffMatrix, LocalStiffMatrix, &
ForceVector, LocalForce, n, NDOFs, DofIndexes, RotateNT, UElement, &
GlobalValues )
!------------------------------------------------------------------------------
TYPE(Matrix_t), POINTER :: StiffMatrix !< The global matrix
REAL(KIND=dp) :: LocalStiffMatrix(:,:) !< Local matrix to be added to the global matrix.
REAL(KIND=dp) :: LocalForce(:) !< Element local force vector.
REAL(KIND=dp) :: ForceVector(:) !< The global RHS vector.
INTEGER :: n !< Number of degrees of freedom in element for each component
INTEGER :: NDOFs !< Number of components for vector field.
INTEGER :: DofIndexes(:) !< Element node/edge/face to global node/edge/face numbering mapping.
LOGICAL, OPTIONAL :: RotateNT !< Should the global equation be done in local normal-tangential coordinates.
TYPE(Element_t), OPTIONAL, TARGET :: UElement !< Element to be updated
REAL(KIND=dp), OPTIONAL :: GlobalValues(:)
!------------------------------------------------------------------------------
INTEGER :: i,j,k,np,dim, NormalIndexes(n), pIndexes(64)
LOGICAL :: Rotate
TYPE(Element_t), POINTER :: Element
TYPE(NormalTangential_t), POINTER :: NT
!------------------------------------------------------------------------------
! Update global matrix and rhs vector....
!------------------------------------------------------------------------------
IF (PRESENT(UElement)) THEN
Element => UElement
ELSE
Element => CurrentModel % CurrentElement
END IF
!------------------------------------------------------------------------------
! Check first if this element has been defined passive
!------------------------------------------------------------------------------
IF ( CheckPassiveElement(Element) ) RETURN
!------------------------------------------------------------------------------
Rotate = .TRUE.
IF ( PRESENT(RotateNT) ) Rotate = RotateNT
dim = CoordinateSystemDimension()
IF( ndofs < dim ) Rotate = .FALSE.
IF( Rotate ) THEN
NT => CurrentModel % Solver % NormalTangential
Rotate = ( NT % NormalTangentialNOFNodes > 0 )
END IF
IF ( Rotate ) THEN
NormalIndexes = 0
np = mGetElementDOFs(pIndexes,Element)
np = MIN(np,n)
NormalIndexes(1:np) = NT % BoundaryReorder(pIndexes(1:np))
CALL RotateMatrix( LocalStiffMatrix, LocalForce, n, dim, NDOFs, &
NormalIndexes, NT % BoundaryNormals, NT % BoundaryTangent1, NT % BoundaryTangent2 )
END IF
!------------------------------------------------------------------------------
IF ( ASSOCIATED( StiffMatrix ) ) THEN
SELECT CASE( StiffMatrix % FORMAT )
CASE( MATRIX_CRS )
CALL CRS_GlueLocalMatrix( StiffMatrix,n,NDOFs, &
DofIndexes, LocalStiffMatrix, GlobalValues )
CASE( MATRIX_LIST )
CALL List_GlueLocalMatrix( StiffMatrix % ListMatrix,n,NDOFs,DofIndexes, &
LocalStiffMatrix )
CASE( MATRIX_BAND,MATRIX_SBAND )
CALL Band_GlueLocalMatrix( StiffMatrix,n,NDOFs,DofIndexes, &
LocalStiffMatrix )
END SELECT
END IF
DO i=1,n
IF ( DofIndexes(i) > 0 ) THEN
DO j=1,NDOFs
k = NDOFs * (DofIndexes(i)-1) + j
!$omp atomic
ForceVector(k) = ForceVector(k) + LocalForce(NDOFs*(i-1)+j)
END DO
END IF
END DO
!------------------------------------------------------------------------------
END SUBROUTINE UpdateGlobalEquations
!------------------------------------------------------------------------------
!> Add element local matrices & vectors to global matrices and vectors.
!> Vectorized version, does not support normal or tangential boundary
!> conditions yet.
SUBROUTINE UpdateGlobalEquationsVec( Gmtr, Lmtr, Gvec, Lvec, n, &
NDOFs, DofIndexes, RotateNT, UElement, MCAssembly )
TYPE(Matrix_t), POINTER :: Gmtr !< The global matrix
REAL(KIND=dp) CONTIG :: Lmtr(:,:) !< Local matrix to be added to the global matrix.
REAL(KIND=dp) CONTIG :: Gvec(:) !< Element local force vector.
REAL(KIND=dp) CONTIG :: Lvec(:) !< The global RHS vector.
INTEGER :: n !< Number of nodes.
INTEGER :: NDOFs !< Number of degrees of free per node.
INTEGER CONTIG :: DofIndexes(:) !< Element node to global node numbering mapping.
LOGICAL, OPTIONAL :: RotateNT !< Should the global equation be done in local normal-tangential coordinates.
TYPE(Element_t), OPTIONAL, TARGET :: UElement !< Element to be updated
LOGICAL, OPTIONAL :: MCAssembly !< Assembly process is multicoloured and guaranteed race condition free
! Local variables
INTEGER :: dim, i,j,k,np
INTEGER :: Ind(n*NDOFs),pIndexes(64)
REAL(KIND=dp) :: Vals(n*NDOFs)
!DIR$ ATTRIBUTES ALIGN:64::Ind, Vals
TYPE(Element_t), POINTER :: Element
LOGICAL :: Rotate
LOGICAL :: ColouredAssembly, NeedMasking
TYPE(NormalTangential_t), POINTER :: NT
IF (PRESENT(UElement)) THEN
Element => UElement
ELSE
Element => CurrentModel % CurrentElement
END IF
IF ( CheckPassiveElement(Element) ) RETURN
Rotate = .TRUE.
IF ( PRESENT(RotateNT) ) Rotate = RotateNT
ColouredAssembly = .FALSE.
IF ( PRESENT(MCAssembly) ) ColouredAssembly = MCAssembly
dim = CoordinateSystemDimension()
IF(ndofs < dim) Rotate = .FALSE.
IF( Rotate ) THEN
NT => CurrentModel % Solver % NormalTangential
Rotate = ( NT % NormalTangentialNOFNodes > 0 )
END IF
IF ( Rotate ) THEN
Ind = 0
np = mGetElementDOFs(pIndexes,Element)
np = MIN(np,n)
Ind(1:np) = NT % BoundaryReorder(pIndexes(1:np))
! TODO: See that RotateMatrix is vectorized
CALL RotateMatrix( Lmtr, Lvec, n, dim, NDOFs, Ind, NT % BoundaryNormals, &
NT % BoundaryTangent1, NT % BoundaryTangent2 )
!IF ( Rotate .AND. NormalTangentialNOFNodes > 0 .AND. ndofs>=dim) THEN
! CALL Fatal('UpdateGlobalEquationsVec', &
! 'Normal or tangential boundary conditions not supported yet!')
END IF
NeedMasking = .FALSE.
DO i=1,n
IF (DofIndexes(i)<=0) THEN
NeedMasking = .TRUE.
EXIT
END IF
END DO
IF ( ASSOCIATED( Gmtr ) ) THEN
SELECT CASE( Gmtr % FORMAT )
CASE( MATRIX_CRS )
CALL CRS_GlueLocalMatrixVec(Gmtr, n, NDOFs, DofIndexes, Lmtr, ColouredAssembly, NeedMasking)
CASE DEFAULT
CALL Fatal('UpdateGlobalEquationsVec','Not implemented for given matrix type')
END SELECT
END IF
! Check for multicolored assembly
IF (ColouredAssembly) THEN
IF (NeedMasking) THEN
! Vector masking needed, no ATOMIC needed
!_ELMER_OMP_SIMD PRIVATE(j,k)
DO i=1,n
IF (DofIndexes(i)>0) THEN
DO j=1,NDOFs
k = NDOFs*(DofIndexes(i)-1) + j
Gvec(k) = Gvec(k) + Lvec(NDOFs*(i-1)+j)
END DO
END IF
END DO
ELSE
! No vector masking needed, no ATOMIC needed
IF (NDOFS>1) THEN
!_ELMER_OMP_SIMD PRIVATE(j,k)
DO i=1,n
DO j=1,NDOFs
k = NDOFs*(DofIndexes(i)-1) + j
Gvec(k) = Gvec(k) + Lvec(NDOFs*(i-1)+j)
END DO
END DO
ELSE
!_ELMER_OMP_SIMD
DO i=1,n
Gvec(DofIndexes(i)) = Gvec(DofIndexes(i)) + Lvec(i)
END DO
END IF
END IF ! Vector masking
ELSE
IF (NeedMasking) THEN
! Vector masking needed, ATOMIC needed
DO i=1,n
IF (DofIndexes(i)>0) THEN
!DIR$ IVDEP
DO j=1,NDOFs
k = NDOFs*(DofIndexes(i)-1) + j
!$OMP ATOMIC
Gvec(k) = Gvec(k) + Lvec(NDOFs*(i-1)+j)
END DO
END IF
END DO
ELSE
! No vector masking needed, ATOMIC needed
DO i=1,n
!DIR$ IVDEP
DO j=1,NDOFs
k = NDOFs*(DofIndexes(i)-1) + j
!$OMP ATOMIC
Gvec(k) = Gvec(k) + Lvec(NDOFs*(i-1)+j)
END DO
END DO
END IF ! Vector masking
END IF ! Coloured assembly
END SUBROUTINE UpdateGlobalEquationsVec
!------------------------------------------------------------------------------
!> Update the global vector with the local vector entry.
!------------------------------------------------------------------------------
SUBROUTINE UpdateGlobalForce(ForceVector, LocalForce, n, &
NDOFs, DofIndexes, RotateNT, UElement )
!------------------------------------------------------------------------------
REAL(KIND=dp) :: LocalForce(:) !< Element local force vector.
REAL(KIND=dp) :: ForceVector(:) !< The global RHS vector.
INTEGER :: n !< Number of nodes.
INTEGER :: NDOFs !< Number of element nodes.
INTEGER :: DofIndexes(:) !< Element node/edge/face to global node/edge/face numbering mapping.
LOGICAL, OPTIONAL :: RotateNT !< Should the global equation be done in local normal-tangential coordinates.
TYPE(Element_t), OPTIONAL, TARGET :: UElement !< Element to be updated
!------------------------------------------------------------------------------
TYPE(Element_t), POINTER :: Element
INTEGER :: i,j,k,np,dim,NormalIndexes(n),pIndexes(64)
LOGICAL :: Rotate
REAL(KIND=dp) :: LocalStiffMatrix(n*NDOFs,n*NDOFs), LForce(n*NDOFs)
TYPE(NormalTangential_t), POINTER :: NT
!------------------------------------------------------------------------------
! Update global matrix and rhs vector....
!------------------------------------------------------------------------------
IF (PRESENT(UElement)) THEN
Element => UElement
ELSE
Element => CurrentModel % CurrentElement
END IF
IF ( CheckPassiveElement( Element ) ) RETURN
Rotate = .TRUE.
IF ( PRESENT(RotateNT) ) Rotate=RotateNT
IF( Rotate ) THEN
NT => CurrentModel % Solver % NormalTangential
Rotate = ( NT % NormalTangentialNOFNodes > 0 )
END IF
IF ( Rotate ) THEN
dim = CoordinateSystemDimension()
NormalIndexes = 0
np = mGetElementDOFs(pIndexes,Element)
np = MIN(np,n)
NormalIndexes(1:np) = NT % BoundaryReorder(pIndexes(1:np))
CALL RotateMatrix( LocalStiffMatrix, LocalForce, n, dim, NDOFs, &
NormalIndexes, NT % BoundaryNormals, NT % BoundaryTangent1, NT % BoundaryTangent2 )
END IF
DO i=1,n
IF ( DofIndexes(i) > 0 ) THEN
DO j=1,NDOFs
k = NDOFs * (DofIndexes(i)-1) + j
!$omp atomic
ForceVector(k) = ForceVector(k) + LocalForce(NDOFs*(i-1)+j)
END DO
END IF
END DO
!------------------------------------------------------------------------------
END SUBROUTINE UpdateGlobalForce
!------------------------------------------------------------------------------