-
Notifications
You must be signed in to change notification settings - Fork 303
/
ElmerSolver.F90
3889 lines (3236 loc) · 155 KB
/
ElmerSolver.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
! *
! *****************************************************************************/
!
!/******************************************************************************
! *
! * ELMER/FEM Solver main program
! *
! ******************************************************************************
! *
! * 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: 02 Jun 1997
! *
! *****************************************************************************/
!> \defgroup Solvers Dynamically linked solvers
!> \defgroup UDF Dynamically linked functions
!> \defgroup Programs Utility programs
!> \degroup ElmerLib Elmer library routines
!> \ingroup ElmerLib
!> \{
#include "../config.h"
!------------------------------------------------------------------------------
!> The main program for Elmer. Solves the equations as defined by the input files.
!------------------------------------------------------------------------------
SUBROUTINE ElmerSolver(initialize)
!------------------------------------------------------------------------------
USE Lists
USE MainUtils
USE OptimizationUtils
USE SolverUtils, ONLY: GetControlValue
!------------------------------------------------------------------------------
IMPLICIT NONE
!------------------------------------------------------------------------------
INTEGER :: Initialize
!------------------------------------------------------------------------------
! Local variables
!------------------------------------------------------------------------------
INTEGER :: i,j,k,n,l,t,k1,k2,iter,Ndeg,istat,nproc,tlen,nthreads
CHARACTER(LEN=MAX_STRING_LEN) :: threads
CHARACTER(:), ALLOCATABLE :: CoordTransform
REAL(KIND=dp) :: s,dt,dtfunc
REAL(KIND=dP), POINTER :: WorkA(:,:,:) => NULL()
REAL(KIND=dp), POINTER, SAVE :: sTime(:), sStep(:), sInterval(:), sSize(:), &
steadyIt(:),nonlinIt(:),sPrevSizes(:,:),sPeriodicTime(:),sPeriodicCycle(:),&
sScan(:),sSweep(:),sPar(:),sFinish(:),sProduce(:),sSlice(:),sSliceRatio(:),&
sSliceWeight(:), sAngle(:), sAngleVelo(:),sSector(:)
LOGICAL :: GotIt,Transient,Scanning, LastSaved, MeshMode = .FALSE.
INTEGER :: TimeIntervals,interval,timestep, &
TotalTimesteps,SavedSteps,CoupledMaxIter,CoupledMinIter
INTEGER, POINTER, SAVE :: Timesteps(:),OutputIntervals(:) => NULL(), ActiveSolvers(:)
REAL(KIND=dp), POINTER, SAVE :: TimestepSizes(:,:),TimestepRatios(:,:)
INTEGER(KIND=AddrInt) :: ControlProcedure
LOGICAL :: InitDirichlet, ExecThis, GotTimestepRatios = .FALSE.
TYPE(ElementType_t),POINTER :: elmt
TYPE(ParEnv_t), POINTER :: ParallelEnv
CHARACTER(LEN=MAX_NAME_LEN) :: ModelName, eq
CHARACTER(LEN=MAX_STRING_LEN) :: OptionString
CHARACTER(:), ALLOCATABLE :: str, PostFile, ExecCommand, OutputFile, RestartFile, &
OutputName, PostName, When
TYPE(Variable_t), POINTER :: Var
TYPE(Mesh_t), POINTER :: Mesh
TYPE(Solver_t), POINTER :: Solver, iSolver
REAL(KIND=dp) :: CT0,RT0,tt
LOGICAL :: Silent=.FALSE., Version=.FALSE., GotModelName, FinishEarly=.FALSE.
LOGICAL :: FirstLoad = .TRUE., FirstTime=.TRUE., Found
INTEGER :: iargc, NoArgs
INTEGER :: iostat, iSweep = 1, OptimIters
LOGICAL :: GotOptimIters
INTEGER :: MeshIndex
TYPE(Mesh_t), POINTER :: ExtrudedMesh
TYPE(Model_t), POINTER, SAVE :: Control
LOGICAL :: DoControl=.FALSE., ProcControl=.FALSE., GotParams=.FALSE., DoIt
INTEGER :: nr,ni,ExtMethod
INTEGER, ALLOCATABLE :: ipar(:)
REAL(KIND=dp), ALLOCATABLE :: rpar(:)
CHARACTER(LEN=MAX_NAME_LEN) :: MeshDir, MeshName
#ifdef HAVE_TRILINOS
INTERFACE
SUBROUTINE TrilinosCleanup() BIND(C,name='TrilinosCleanup')
IMPLICIT NONE
END SUBROUTINE TrilinosCleanup
END INTERFACE
#endif
! Start the watches, store later
!--------------------------------
RT0 = RealTime()
CT0 = CPUTime()
! If parallel execution requested, initialize parallel environment:
!------------------------------------------------------------------
IF(FirstTime) ParallelEnv => ParallelInit()
OutputPE = -1
IF( ParEnv % MyPe == 0 ) THEN
OutputPE = 0
END IF
IF ( FirstTime ) THEN
!
! Print banner to output:
! -----------------------
NoArgs = COMMAND_ARGUMENT_COUNT()
! Info Level is always true until the model has been read!
! This makes it possible to cast something
Silent = .FALSE.
Version = .FALSE.
IF( NoArgs > 0 ) THEN
i = 0
DO WHILE( i < NoArgs )
i = i + 1
CALL GET_COMMAND_ARGUMENT(i, OptionString)
IF( OptionString=='-rpar' ) THEN
! Followed by number of parameters + the parameter values
i = i + 1
CALL GET_COMMAND_ARGUMENT(i, OptionString)
READ( OptionString,*) nr
ALLOCATE( rpar(nr) )
DO j=1,nr
i = i + 1
CALL GET_COMMAND_ARGUMENT(i, OptionString)
READ( OptionString,*) rpar(j)
END DO
CALL Info('MAIN','Read '//I2S(nr)//' real parameters from command line!')
CALL SetRealParametersMATC(nr,rpar)
END IF
IF( OptionString=='-ipar' ) THEN
! Followed by number of parameters + the parameter values
i = i + 1
CALL GET_COMMAND_ARGUMENT(i, OptionString)
READ( OptionString,*) ni
ALLOCATE( ipar(nr) )
DO j=1,ni
i = i + 1
CALL GET_COMMAND_ARGUMENT(i, OptionString)
READ( OptionString,*) ipar(j)
END DO
CALL Info('MAIN','Read '//I2S(ni)//' integer parameters from command line!')
CALL SetIntegerParametersMATC(ni,ipar)
END IF
Silent = Silent .OR. &
( OptionString=='-s' .OR. OptionString=='--silent' )
Version = Version .OR. &
( OptionString=='-v' .OR. OptionString == '--version' )
END DO
END IF
! Set number of OpenMP threads
nthreads = 1
!$ nthreads = omp_get_max_threads()
IF (nthreads > 1 ) THEN
! Check if OMP_NUM_THREADS environment variable is set
CALL envir( 'OMP_NUM_THREADS', threads, tlen )
IF (tlen==0) THEN
CALL Info('MAIN','OMP_NUM_THREADS not set. Using only 1 thread per task.',Level=6)
nthreads = 1
! Set number of threads to 1
!$ CALL omp_set_num_threads(nthreads)
#ifdef HAVE_MKL
CALL mkl_set_num_threads(nthreads)
#endif
END IF
END IF
ParEnv % NumberOfThreads = nthreads
IF( .NOT. Silent ) THEN
CALL Info( 'MAIN', ' ')
CALL Info( 'MAIN', '=============================================================')
CALL Info( 'MAIN', 'ElmerSolver finite element software, Welcome! ')
CALL Info( 'MAIN', 'This program is free software licensed under (L)GPL ')
CALL Info( 'MAIN', 'Copyright 1st April 1995 - , CSC - IT Center for Science Ltd.')
CALL Info( 'MAIN', 'Webpage http://www.csc.fi/elmer, Email elmeradm@csc.fi ')
CALL Info( 'MAIN', 'Version: ' // GetVersion() // ' (Rev: ' // GetRevision() // &
', Compiled: ' // GetCompilationDate() // ')' )
IF ( ParEnv % PEs > 1 ) THEN
CALL Info( 'MAIN', ' Running in parallel using ' // &
i2s(ParEnv % PEs) // ' tasks.')
ELSE
CALL Info('MAIN', ' Running one task without MPI parallelization.',Level=10)
END IF
! Print out number of threads in use
IF ( nthreads > 1 ) THEN
CALL Info('MAIN', ' Running in parallel with ' // &
i2s(nthreads) // ' threads per task.')
ELSE
CALL Info('MAIN', ' Running with just one thread per task.',Level=10)
END IF
#ifdef HAVE_FETI4I
CALL Info( 'MAIN', ' FETI4I library linked in.')
#endif
#ifdef HAVE_HYPRE
CALL Info( 'MAIN', ' HYPRE library linked in.')
#endif
#ifdef HAVE_TRILINOS
CALL Info( 'MAIN', ' Trilinos library linked in.')
#endif
#ifdef HAVE_MUMPS
CALL Info( 'MAIN', ' MUMPS library linked in.')
#endif
#ifdef HAVE_CHOLMOD
CALL Info( 'MAIN', ' CHOLMOD library linked in.')
#endif
#ifdef HAVE_SUPERLU
CALL Info( 'MAIN', ' SUPERLU library linked in.')
#endif
#ifdef HAVE_PARDISO
CALL Info( 'MAIN', ' PARDISO library linked in.')
#endif
#ifdef HAVE_MMG
CALL Info( 'MAIN', ' MMG library linked in.')
#endif
#ifdef HAVE_PARMMG
CALL Info( 'MAIN', ' ParMMG library linked in.')
#endif
#ifdef HAVE_MKL
CALL Info( 'MAIN', ' Intel MKL linked in.' )
#endif
#ifdef HAVE_LUA
CALL Info( 'MAIN', ' Lua interpreter linked in.' )
#endif
#ifdef HAVE_EXTOPTIM
CALL Info( 'MAIN', ' External optimization routines linked in.' )
#endif
#ifdef HAVE_ZOLTAN
CALL Info( 'MAIN', ' Zoltan library linked in.' )
#endif
#ifdef HAVE_AMGX
CALL Info( 'MAIN', ' AMGX library linked in.' )
#endif
CALL Info( 'MAIN', '=============================================================')
END IF
IF( Version ) RETURN
CALL InitializeElementDescriptions()
END IF
! Read input file name either as an argument, or from the default file:
!----------------------------------------------------------------------
GotModelName = .FALSE.
IF ( NoArgs > 0 ) THEN
CALL GET_COMMAND_ARGUMENT(1, ModelName)
IF( ModelName(1:1) /= '-') THEN
GotModelName = .TRUE.
IF (NoArgs > 1) CALL GET_COMMAND_ARGUMENT(2, eq)
END IF
END IF
IF( .NOT. GotModelName ) THEN
OPEN( 1, File='ELMERSOLVER_STARTINFO', STATUS='OLD', IOSTAT=iostat )
IF( iostat /= 0 ) THEN
CALL Fatal( 'MAIN', 'Unable to find ELMERSOLVER_STARTINFO, can not execute.' )
END IF
READ(1,'(a)') ModelName
CLOSE(1)
END IF
! This sets optionally some internal parameters for doing scanning
! over a parameter space / optimization.
!-----------------------------------------------------------------
IF( FirstTime ) THEN
OPEN( Unit=InFileUnit, Action='Read',File=ModelName,Status='OLD',IOSTAT=iostat)
IF( iostat /= 0 ) THEN
CALL Fatal( 'MAIN', 'Unable to find input file [' // &
TRIM(Modelname) // '], can not execute.' )
END IF
ALLOCATE( Control )
! Read only the "Run Control" section of the sif file.
CALL LoadInputFile( Control,InFileUnit,ModelName,MeshDir,MeshName, &
.FALSE., .TRUE., ControlOnly = .TRUE.)
DoControl = ASSOCIATED( Control % Control )
IF( DoControl ) THEN
CALL Info('MAIN','Run Control section active!')
OptimIters = ListGetInteger( Control % Control,'Run Control Iterations', GotOptimIters )
IF(.NOT. GotOptimIters) OptimIters = 1
! If there are no parameters this does nothing
CALL ControlParameters(Control % Control,1,GotParams,FinishEarly)
ELSE
OptimIters = 1
END IF
END IF
!------------------------------------------------------------------------------
! Read element definition file, and initialize element types
!------------------------------------------------------------------------------
IF ( Initialize==1 ) THEN
CALL FreeModel(CurrentModel)
FirstLoad=.TRUE.
END IF
!------------------------------------------------------------------------------
! Read Model and mesh from Elmer mesh data base
!------------------------------------------------------------------------------
MeshIndex = 0
DO WHILE( .TRUE. )
IF ( initialize==2 ) GOTO 1
IF(MeshMode) THEN
CALL FreeModel(CurrentModel)
MeshIndex = MeshIndex + 1
FirstLoad = .TRUE.
END IF
IF ( FirstLoad ) THEN
IF( .NOT. Silent ) THEN
CALL Info( 'MAIN', ' ')
CALL Info( 'MAIN', ' ')
CALL Info( 'MAIN', '-------------------------------------')
CALL Info( 'MAIN', 'Reading Model: '//TRIM( ModelName) )
END IF
INQUIRE(Unit=InFileUnit, Opened=GotIt)
IF ( gotIt ) CLOSE(inFileUnit)
! Here we read the whole model including command file and default mesh file
!---------------------------------------------------------------------------------
OPEN( Unit=InFileUnit, Action='Read',File=ModelName,Status='OLD',IOSTAT=iostat)
IF( iostat /= 0 ) THEN
CALL Fatal( 'MAIN', 'Unable to find input file [' // &
TRIM(Modelname) // '], can not execute.' )
END IF
CurrentModel => LoadModel(ModelName,.FALSE.,ParEnv % PEs,ParEnv % MyPE,MeshIndex)
IF(.NOT.ASSOCIATED(CurrentModel)) EXIT
!----------------------------------------------------------------------------------
! Set namespace searching mode
!----------------------------------------------------------------------------------
CALL SetNamespaceCheck( ListGetLogical( CurrentModel % Simulation, &
'Additive namespaces', Found ) )
!----------------------------------------------------------------------------------
MeshMode = ListGetLogical( CurrentModel % Simulation, 'Mesh Mode', Found)
!------------------------------------------------------------------------------
! Some keywords automatically require other keywords to be set
! We could complain on the missing keywords later on, but sometimes
! it may be just as simple to add them directly.
!------------------------------------------------------------------------------
CALL CompleteModelKeywords( )
!----------------------------------------------------------------------------------
! Optionally perform simple extrusion to increase the dimension of the mesh
!----------------------------------------------------------------------------------
CALL CreateExtrudedMesh()
!----------------------------------------------------------------------------------
! If requested perform coordinate transformation directly after is has been obtained.
! Don't maintain the original mesh.
!----------------------------------------------------------------------------------
CoordTransform=ListGetString(CurrentModel % Simulation,'Coordinate Transformation',GotIt)
IF( GotIt ) THEN
CALL CoordinateTransformation( CurrentModel % Meshes, CoordTransform, &
CurrentModel % Simulation, .TRUE. )
END IF
IF(.NOT. Silent ) THEN
CALL Info( 'MAIN', '-------------------------------------')
END IF
ELSE
IF ( Initialize==3 ) THEN
INQUIRE(Unit=InFileUnit, Opened=GotIt)
IF ( gotIt ) CLOSE(inFileUnit)
OPEN( Unit=InFileUnit, Action='Read', &
File=ModelName,Status='OLD',IOSTAT=iostat)
IF( iostat /= 0 ) THEN
CALL Fatal( 'MAIN', 'Unable to find input file [' // &
TRIM(Modelname) // '], can not execute.' )
END IF
END IF
IF ( .NOT.ReloadInputFile(CurrentModel) ) EXIT
Mesh => CurrentModel % Meshes
DO WHILE( ASSOCIATED(Mesh) )
Mesh % SavesDone = 0
Mesh => Mesh % Next
END DO
END IF
1 CONTINUE
CALL ListAddLogical( CurrentModel % Simulation, &
'Initialization Phase', .TRUE. )
! Save the start time to the list so that it may be retrieved when necessary
! This could perhaps also be a global variable etc, but this will do for now.
!-------------------------------------------------------------------------
IF( ListGetLogical( CurrentModel % Simulation,'Simulation Timing',GotIt) ) THEN
CALL ListAddConstReal( CurrentModel % Simulation,'cputime0',ct0 )
CALL ListAddConstReal( CurrentModel % Simulation,'realtime0',rt0 )
END IF
!------------------------------------------------------------------------------
! Check for transient case
!------------------------------------------------------------------------------
eq = ListGetString( CurrentModel % Simulation, 'Simulation Type', GotIt )
Scanning = ( eq == 'scanning' )
Transient = ( eq == 'transient' )
!------------------------------------------------------------------------------
! To more conveniently support the use of VTK based visualization there
! is a hack that recognizes VTU suffix and creates a instance of output solver.
! Note that this is really quite a dirty hack, and is not a good example.
!-----------------------------------------------------------------------------
CALL AddVtuOutputSolverHack()
! Support easily saving scalars to a file activated by "Scalars File" keyword.
!-----------------------------------------------------------------------------
CALL AddSaveScalarsHack()
!------------------------------------------------------------------------------
! Add coordinates such that if there is a solver that is run on creation
! the coordinates are already usable then.
!------------------------------------------------------------------------------
IF ( FirstLoad ) CALL AddMeshCoordinates()
!------------------------------------------------------------------------------
! Figure out what (flow,heat,stress,...) should be computed, and get
! memory for the dofs
!------------------------------------------------------------------------------
CALL AddSolvers()
!------------------------------------------------------------------------------
! Time integration and/or steady state steps
!------------------------------------------------------------------------------
CALL InitializeIntervals()
!------------------------------------------------------------------------------
! Add time and other global variables so that we can have dependence on these.
!------------------------------------------------------------------------------
IF ( FirstLoad ) CALL AddTimeEtc()
!------------------------------------------------------------------------------
! Initialize the random seeds so that all simulation depending on that
! give consistent results.
!------------------------------------------------------------------------------
IF( FirstLoad ) CALL InitializeRandomSeed()
!------------------------------------------------------------------------------
! Get Output File Options
!------------------------------------------------------------------------------
! Initial Conditions:
! -------------------
IF ( FirstLoad ) THEN
CALL SetInitialConditions()
DO i=1,CurrentModel % NumberOfSolvers
Solver => CurrentModel % Solvers(i)
IF( ListGetLogical( Solver % Values, 'Initialize Exported Variables', GotIt ) ) THEN
CurrentModel % Solver => Solver
CALL UpdateExportedVariables( Solver )
END IF
END DO
END IF
! Compute the total number of steps that will be saved to the files
! Particularly look if the last step will be saved, or if it has
! to be saved separately.
!------------------------------------------------------------------
CALL CountSavedTimesteps()
CALL ListAddLogical( CurrentModel % Simulation, &
'Initialization Phase', .FALSE. )
FirstLoad = .FALSE.
IF ( Initialize == 1 ) EXIT
! Check whether we are using external optimization routine that
! needs to have basically Elmer given as a function that returns
! the cost function. Hence this is treated separately from the internal
! optimization methods.
!---------------------------------------------------------------------
ExtMethod = 0
str = ListGetString( CurrentModel % Control,'Optimization Method',Found)
IF( Found ) THEN
IF( SEQL(str, 'hybrd') ) THEN
ExtMethod = 1
ELSE IF( SEQL(str,'newuoa') ) THEN
ExtMethod = 2
ELSE IF( SEQL(str,'bobyqa') ) THEN
ExtMethod = 3
END IF
END IF
ExecCommand = ListGetString( CurrentModel % Simulation, &
'Control Procedure', ProcControl )
IF ( ProcControl ) THEN
ControlProcedure = GetProcAddr( ExecCommand )
CALL ExecSimulationProc( ControlProcedure, CurrentModel )
ELSE IF( ExtMethod > 0 ) THEN
#ifdef HAVE_EXTOPTIM
SELECT CASE( ExtMethod )
CASE(1)
CALL ExternalOptimization_minpack(ExecSimulationFunVec)
CASE(2)
CALL ExternalOptimization_newuoa(ExecSimulationFunCost)
CASE(3)
CALL ExternalOptimization_bobyqa(ExecSimulationFunCost)
END SELECT
#else
CALL Fatal('MAIN','Compile WITH_EXTOPTIM to activate method: '//TRIM(str))
#endif
ELSE IF( DoControl ) THEN
! This sets optionally some internal parameters for doing scanning
! over a parameter space / optimization.
!-----------------------------------------------------------------
iSweep = 0
DO WHILE (.TRUE.)
iSweep = iSweep + 1
CALL Info('MAIN','========================================================',Level=5)
CALL Info('MAIN','Control Loop '//I2S(iSweep))
CALL Info('MAIN','========================================================',Level=5)
sSweep = 1.0_dp * iSweep
! If there are no parameters this does nothing
CALL ControlResetMesh(CurrentModel % Control, iSweep )
IF( iSweep > 1 ) THEN
CALL ControlParameters(CurrentModel % Control,iSweep,&
GotParams,FinishEarly)
IF( FinishEarly ) EXIT
Found = ReloadInputFile(CurrentModel,RewindFile=.TRUE.)
CALL InitializeIntervals()
END IF
! This is another calling slot as here we have formed the model structure and
! may toggle with the keyword coefficients.
CALL ControlParameters(CurrentModel % Control,iSweep,&
GotParams,FinishEarly,SetCoeffs=.TRUE.)
IF( iSweep > 1 ) THEN
IF( ListGetLogical( CurrentModel % Control,'Reset Adaptive Mesh',Found ) ) THEN
CALL ResetAdaptiveMesh()
END IF
IF( ListGetLogical( CurrentModel % Control,'Reset Initial Conditions',Found ) ) THEN
CALL SetInitialConditions()
END IF
END IF
CALL ExecSimulation( TimeIntervals, CoupledMinIter, &
CoupledMaxIter, OutputIntervals, Transient, Scanning)
! This evaluates the cost function and saves the results of control
CALL ControlParameters(CurrentModel % Control, &
iSweep,GotParams,FinishEarly,.TRUE.)
IF( iSweep == 1 ) THEN
DO i=1,CurrentModel % NumberOfSolvers
iSolver => CurrentModel % Solvers(i)
j = iSolver % NumberOfConstraintModes
IF( j <= 0 ) CYCLE
IF( ListGetLogical( iSolver % Values,'Run Control Constraint Modes', Found ) .OR. &
ListGetLogical( CurrentModel % Control,'Constraint Modes Analysis',Found ) ) THEN
IF( GotOptimIters ) THEN
IF( OptimIters /= j ) THEN
CALL Warn('MAIN','Incompatible number of run control iterations and constraint modes!')
END IF
ELSE
CALL Info('MAIN','Setting run control iterations to constraint modes count: '//I2S(j))
OptimIters = j
END IF
EXIT
END IF
END DO
END IF
! We use this type of condition so that OptimIters can be changed on-the-fly
IF(iSweep == OptimIters) EXIT
END DO
DO i=1,CurrentModel % NumberOfSolvers
iSolver => CurrentModel % Solvers(i)
IF( iSolver % NumberOfConstraintModes > 0 ) THEN
IF( ListGetLogical( iSolver % Values,'Run Control Constraint Modes', Found ) .OR. &
ListGetLogical( CurrentModel % Control,'Constraint Modes Analysis',Found ) ) THEN
CALL FinalizeLumpedMatrix( iSolver )
END IF
END IF
END DO
DO i=1,CurrentModel % NumberOfSolvers
iSolver => CurrentModel % Solvers(i)
IF ( iSolver % PROCEDURE == 0 ) CYCLE
When = ListGetString( iSolver % Values, 'Exec Solver', Found )
IF ( Found ) THEN
DoIt = ( When == 'after control' )
ELSE
DoIt = ( iSolver % SolverExecWhen == SOLVER_EXEC_AFTER_CONTROL )
END IF
IF(DoIt) CALL SolverActivate( CurrentModel,iSolver,dt,Transient )
END DO
ELSE
CALL ExecSimulation( TimeIntervals, CoupledMinIter, &
CoupledMaxIter, OutputIntervals, Transient, Scanning)
END IF
! Comparison to reference is done to enable consistency test under ctest.
!-------------------------------------------------------------------------
CALL CompareToReferenceSolution( )
IF ( Initialize >= 2 ) EXIT
END DO
IF( ListGetLogical( CurrentModel % Simulation,'Echo Keywords at End', GotIt ) ) THEN
CALL ListEchoKeywords( CurrentModel )
END IF
CALL CompareToReferenceSolution( Finalize = .TRUE. )
#ifdef DEVEL_LISTUSAGE
IF(InfoActive(10) .AND. ParEnv % MyPe == 0 ) THEN
CALL Info('MAIN','Reporting unused list entries for sif improvement!')
CALL Info('MAIN','If you do not want these lines undefine > DEVEL_LISTUSAGE < !')
CALL ReportListCounters( CurrentModel, 2 )
END IF
#endif
#ifdef DEVEL_LISTCOUNTER
IF( ParEnv % MyPe == 0 ) THEN
CALL Info('MAIN','Reporting list counters for code optimization purposes only!')
CALL Info('MAIN','If you get these lines with production code undefine > DEVEL_LISTCOUNTER < !')
CALL ReportListCounters( CurrentModel, 3 )
END IF
#endif
!------------------------------------------------------------------------------
! THIS IS THE END (...,at last, the end, my friend,...)
!------------------------------------------------------------------------------
IF ( Initialize /= 1 ) CALL Info( 'MAIN', '*** Elmer Solver: ALL DONE ***',Level=3 )
! This may be used to study problems at the finish
IF( ListGetLogical( CurrentModel % Simulation,'Dirty Finish', GotIt ) ) THEN
CALL Info('MAIN','Skipping freeing of the Model structure',Level=4)
RETURN
END IF
IF ( Initialize <= 0 ) CALL FreeModel(CurrentModel)
#ifdef HAVE_TRILINOS
CALL TrilinosCleanup()
#endif
IF ( FirstTime ) CALL ParallelFinalize()
FirstTime = .FALSE.
CALL Info('MAIN','The end',Level=3)
RETURN
CONTAINS
! If we want to start a new adaptive simulation with the original mesh
! call this subroutine.
!---------------------------------------------------------------------
SUBROUTINE ResetAdaptiveMesh()
TYPE(Mesh_t), POINTER :: pMesh, pMesh0
TYPE(Solver_t), POINTER :: iSolver
TYPE(Variable_t), POINTER :: pVar
LOGICAL :: GB, BO
CHARACTER(*), PARAMETER :: Caller = 'ResetAdaptiveMesh'
! Find the 1st mesh
pMesh0 => CurrentModel % Mesh
DO WHILE( ASSOCIATED(pMesh0 % Parent) )
pMesh0 => pMesh0 % Parent
END DO
!PRINT *,'First mesh:',pMesh0 % AdaptiveDepth, TRIM(pMesh0 % Name)
! Find the last mesh
pMesh => CurrentModel % Mesh
DO WHILE( ASSOCIATED(pMesh % Child) )
pMesh => pMesh % Child
END DO
!PRINT *,'Last mesh:',pMesh % AdaptiveDepth, TRIM(pMesh % Name)
! Move point to the 1st mesh and related fields
CALL SetCurrentMesh( CurrentModel, pMesh0 )
DO i=1,CurrentModel % NumberOfSolvers
iSolver => CurrentModel % Solvers(i)
IF(.NOT. ASSOCIATED(iSolver % Variable)) THEN
CALL Info(Caller,'No Variable in this mesh for solver index '//I2S(i),Level=10)
CYCLE
END IF
! Set Solver mesh
IF(ASSOCIATED(iSolver % Mesh)) iSolver % Mesh => pMesh0
! Set Solver variable point to the field in the original mesh
IF(ASSOCIATED(iSolver % Variable)) THEN
pVar => VariableGet(pMesh0 % Variables, &
iSolver % Variable % Name, ThisOnly = .TRUE.)
IF(.NOT. ASSOCIATED(pVar)) THEN
CALL Info(Caller,'No Variable in coarsest mesh for solver index '//I2S(i),Level=10)
CYCLE
END IF
iSolver % Variable => pVar
END IF
! Reset active element table
IF( iSolver % NumberOfActiveElements == 0 ) THEN
CALL Info(Caller,'No active elements for solver index '//I2S(i),Level=10)
CYCLE
END IF
iSolver % NumberOfActiveElements = 0
CALL SetActiveElementsTable( CurrentModel, iSolver )
! Create the matrix related to the original mesh
IF( .NOT. ASSOCIATED( iSolver % Matrix ) ) THEN
CALL Info(Caller,'No matrix for solver index '//I2S(i),Level=10)
CYCLE
END IF
CALL FreeMatrix( iSolver % Matrix)
GB = ListGetLogical( iSolver % Values,'Bubbles in Global System', Found )
IF ( .NOT. Found ) GB = .TRUE.
BO = ListGetLogical( iSolver % Values,'Optimize Bandwidth', Found )
IF ( .NOT. Found ) BO = .TRUE.
iSolver % Matrix => CreateMatrix( CurrentModel, iSolver, iSolver % Mesh, &
iSolver % Variable % Perm, iSolver % Variable % DOFs, MATRIX_CRS, &
BO, ListGetString( iSolver % Values, 'Equation' ), GlobalBubbles=GB )
ALLOCATE( iSolver % Matrix % rhs(iSolver % Matrix % NumberOfRows ) )
iSolver % Matrix % rhs = 0.0_dp
END DO
! Release the old adaptive meshes
DO WHILE( ASSOCIATED(pMesh % Parent))
pMesh => pMesh % Parent
CALL ReleaseMesh( pMesh % Child )
END DO
pMesh % Child => NULL()
END SUBROUTINE ResetAdaptiveMesh
SUBROUTINE InitializeRandomSeed()
INTEGER :: i,n
INTEGER, ALLOCATABLE :: seeds(:)
CALL RANDOM_SEED() ! initialize with system generated seed
i = ListGetInteger( CurrentModel % Simulation,'Random Number Seed',Found )
IF( .NOT. Found ) i = 314159265
CALL RANDOM_SEED(size=j) ! find out size of seed
ALLOCATE(seeds(j))
!CALL RANDOM_SEED(get=seeds) ! get system generated seed
!WRITE(*,*) seeds ! writes system generated seed
seeds = i
CALL RANDOM_SEED(put=seeds) ! set current seed
!CALL RANDOM_SEED(get=seeds) ! get current seed
!WRITE(*,*) seeds ! writes 314159265
DEALLOCATE(seeds)
CALL Info('MAIN','Random seed initialized to: '//I2S(i),Level=10)
END SUBROUTINE InitializeRandomSeed
! Optionally create extruded mesh on-the-fly.
!--------------------------------------------------------------------
SUBROUTINE CreateExtrudedMesh()
INTEGER :: ExtrudeLayers
LOGICAL :: SliceVersion
IF(.NOT. ListCheckPrefix(CurrentModel % Simulation,'Extruded Mesh') ) RETURN
ExtrudeLayers = GetInteger(CurrentModel % Simulation,'Extruded Mesh Levels',Found)-1
IF( .NOT. Found ) THEN
ExtrudeLayers = GetInteger(CurrentModel % Simulation,'Extruded Mesh Layers',Found)
END IF
IF(.NOT. Found ) RETURN
IF(ExtrudeLayers < 2) THEN
CALL Fatal('MAIN','There must be at least two "Extruded Mesh Layers"!')
END IF
SliceVersion = GetLogical(CurrentModel % Simulation,'Extruded Mesh Slices',Found )
IF( SliceVersion ) THEN
ExtrudedMesh => MeshExtrudeSlices(CurrentModel % Meshes, ExtrudeLayers-1)
ELSE
ExtrudedMesh => MeshExtrude(CurrentModel % Meshes, ExtrudeLayers-1)
END IF
! Make the solvers point to the extruded mesh, not the original mesh
!-------------------------------------------------------------------
DO i=1,CurrentModel % NumberOfSolvers
IF(ASSOCIATED(CurrentModel % Solvers(i) % Mesh,CurrentModel % Meshes)) &
CurrentModel % Solvers(i) % Mesh => ExtrudedMesh
END DO
ExtrudedMesh % Next => CurrentModel % Meshes % Next
CurrentModel % Meshes => ExtrudedMesh
! If periodic BC given, compute boundary mesh projector:
! ------------------------------------------------------
DO i = 1,CurrentModel % NumberOfBCs
IF(ASSOCIATED(CurrentModel % Bcs(i) % PMatrix)) &
CALL FreeMatrix( CurrentModel % BCs(i) % PMatrix )
CurrentModel % BCs(i) % PMatrix => NULL()
k = ListGetInteger( CurrentModel % BCs(i) % Values, 'Periodic BC', GotIt )
IF( GotIt ) THEN
CurrentModel % BCs(i) % PMatrix => PeriodicProjector( CurrentModel, ExtrudedMesh, i, k )
END IF
END DO
END SUBROUTINE CreateExtrudedMesh
! Given geometric ratio of timesteps redistribute them so that the ratio
! is met as closely as possible, maintaining total time and sacrificing
! number of timesteps.
!----------------------------------------------------------------------
SUBROUTINE GeometricTimesteps(m,n0,dt0,r)
INTEGER :: m
INTEGER :: n0(:)
REAL(KIND=dp) :: dt0(:),r(:)
INTEGER :: i,n
REAL(KIND=dp) :: q
LOGICAL :: Visited = .FALSE.
! Only do this once since it tampers stuff in lists.
IF(Visited) RETURN
Visited = .TRUE.
CALL Info('MAIN','Creating geometric timestepping strategy',Level=6)
DO i=1,m
! Some users may give zero ratio, assume that they mean one.
IF(ABS(r(i)) < EPSILON(q) ) r(i) = 1.0_dp
! Ratio one means even distribution.
IF(ABS(r(i)-1.0) < EPSILON(q) ) CYCLE
q = 1 + (r(i)-1)/n0(i)
n = NINT( LOG(1+(q-1)*n0(i)) / LOG(q) )
dt = n0(i)*dt0(i)*(1-q)/(1-q**n)
!PRINT *,'ratio:',i,n0(i),dt0(i),r(i),n,dt,q
! Replace the new distribution
r(i) = q
dt0(i) = dt
n0(i) = n
END DO
END SUBROUTINE GeometricTimesteps
! Initialize intervals for steady state, transient and scanning types.
!----------------------------------------------------------------------
SUBROUTINE InitializeIntervals()
IF ( Transient .OR. Scanning ) THEN
Timesteps => ListGetIntegerArray( CurrentModel % Simulation, &
'Timestep Intervals', GotIt )
IF ( .NOT.GotIt ) THEN
CALL Fatal('MAIN', 'Keyword > Timestep Intervals < MUST be ' // &
'defined for transient and scanning simulations' )
END IF
TimestepSizes => ListGetConstRealArray( CurrentModel % Simulation, &
'Timestep Sizes', GotIt )
IF ( .NOT.GotIt ) THEN
IF( Scanning .OR. ListCheckPresent( CurrentModel % Simulation,'Timestep Size') ) THEN
ALLOCATE(TimestepSizes(SIZE(Timesteps),1))
TimestepSizes = 1.0_dp
ELSE
CALL Fatal( 'MAIN', 'Keyword [Timestep Sizes] MUST be ' // &
'defined for time dependent simulations' )
END IF
END IF
CoupledMaxIter = ListGetInteger( CurrentModel % Simulation, &
'Steady State Max Iterations', GotIt, minv=1 )
IF ( .NOT. GotIt ) CoupledMaxIter = 1
TimeIntervals = SIZE(Timesteps)
TimestepRatios => ListGetConstRealArray( CurrentModel % Simulation, &
'Timestep Ratios', GotTimestepRatios )
IF ( GotTimestepRatios ) THEN
CALL GeometricTimesteps(TimeIntervals,Timesteps,TimestepSizes(:,1),TimestepRatios(:,1))
END IF
ELSE
! Steady state
!------------------------------------------------------------------------------
IF( .NOT. ASSOCIATED(Timesteps) ) THEN
ALLOCATE(Timesteps(1))
END IF
Timesteps(1) = ListGetInteger( CurrentModel % Simulation, &
'Steady State Max Iterations', GotIt,minv=1 )
IF ( .NOT. GotIt ) Timesteps(1) = 1
CoupledMaxIter = 1
ALLOCATE(TimestepSizes(1,1))
TimestepSizes(1,1) = 1.0_dp
TimeIntervals = 1
END IF
CoupledMinIter = ListGetInteger( CurrentModel % Simulation, &
'Steady State Min Iterations', GotIt )
IF( .NOT. GotIt ) CoupledMinIter = 1
OutputIntervals => ListGetIntegerArray( CurrentModel % Simulation, &
'Output Intervals', GotIt )
IF( GotIt ) THEN
IF( SIZE(OutputIntervals) /= SIZE(TimeSteps) ) THEN
CALL Fatal('MAIN','> Output Intervals < should have the same size as > Timestep Intervals < !')
END IF
ELSE
IF( .NOT. ASSOCIATED( OutputIntervals ) ) THEN
ALLOCATE( OutputIntervals(SIZE(TimeSteps)) )
OutputIntervals = 1
END IF
END IF