ElmerFoamFSI  2.0
ElmerFoamFSI is fluid-solid interaction simulation application built up from OpenFOAM CFD and Elmer CSM coupled through the IMPACT multiphysics software integration infrastructure.
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Groups Pages
Initialize.F90
Go to the documentation of this file.
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 ! * ELMER/FEM Solver main program
27 ! *
28 ! ******************************************************************************
29 ! *
30 ! * Authors: Juha Ruokolainen
31 ! * Email: Juha.Ruokolainen@csc.fi
32 ! * Web: http://www.csc.fi/elmer
33 ! * Address: CSC - IT Center for Science Ltd.
34 ! * Keilaranta 14
35 ! * 02101 Espoo, Finland
36 ! *
37 ! * Original Date: 02 Jun 1997
38 ! *
39 ! *****************************************************************************/
40 
42 
44 
46 
48 
51 
52 !------------------------------------------------------------------------------
54 !------------------------------------------------------------------------------
55  SUBROUTINE elmerinitialize(runs)
56 !------------------------------------------------------------------------------
57  USE testobject
58  USE generalmodule
59 !------------------------------------------------------------------------------
60  IMPLICIT NONE
61 !------------------------------------------------------------------------------
62 
63 ! INTEGER :: Initialize
64 
65 !------------------------------------------------------------------------------
66 ! Local variables
67 !------------------------------------------------------------------------------
68 
69 ! INTEGER :: i,j,k,n,l,t,k1,k2,iter,Ndeg,istat,nproc,tlen,nthreads
70  CHARACTER(LEN=MAX_STRING_LEN) :: coordtransform
71 ! CHARACTER(LEN=MAX_STRING_LEN) :: threads
72 !
73 ! REAL(KIND=dp) :: s,dt,dtfunc
74 ! REAL(KIND=dP), POINTER :: WorkA(:,:,:) => NULL()
75 ! REAL(KIND=dp), POINTER, SAVE :: sTime(:), sStep(:), sInterval(:), sSize(:), &
76 ! steadyIt(:),nonlinIt(:),sPrevSizes(:,:),sPeriodic(:)
77 !
78 ! TYPE(Element_t),POINTER :: CurrentElement
79 !
80 ! LOGICAL :: GotIt,Transient,Scanning,LastSaved
81 !
82 ! INTEGER :: TimeIntervals,interval,timestep, &
83 ! TotalTimesteps,SavedSteps,CoupledMaxIter,CoupledMinIter
84 !
85 ! INTEGER, POINTER, SAVE :: Timesteps(:),OutputIntervals(:), ActiveSolvers(:)
86 ! REAL(KIND=dp), POINTER, SAVE :: TimestepSizes(:,:)
87 !
88 ! INTEGER(KIND=AddrInt) :: ControlProcedure
89 !
90 ! LOGICAL :: InitDirichlet, ExecThis
91 !
92 ! TYPE(ElementType_t),POINTER :: elmt
93 !
94 ! TYPE(ParEnv_t), POINTER :: ParallelEnv
95 !
96 ! CHARACTER(LEN=MAX_NAME_LEN) :: ModelName, eq, ExecCommand
97 ! CHARACTER(LEN=MAX_STRING_LEN) :: OutputFile, PostFile, RestartFile, &
98 ! OutputName=' ',PostName=' ', When, OptionString
99 !
100 ! TYPE(Variable_t), POINTER :: Var
101 ! TYPE(Mesh_t), POINTER :: Mesh
102 ! TYPE(Solver_t), POINTER :: Solver
103 !
104 ! REAL(KIND=dp) :: RealTime,CPUTime
105 !
106 ! LOGICAL :: FirstLoad = .TRUE., FirstTime=.TRUE., Found
107 ! LOGICAL :: Silent, Version, GotModelName
108 !
109  INTEGER :: iargc
110 !
111 ! INTEGER :: ExtrudeLevels
112 ! TYPE(Mesh_t), POINTER :: ExtrudedMesh
113 !
114 ! INTEGER :: omp_get_max_threads
115 !
116 !#ifdef HAVE_TRILINOS
117 !INTERFACE
118 ! SUBROUTINE TrilinosCleanup() BIND(C,name='TrilinosCleanup')
119 ! IMPLICIT NONE
120 ! END SUBROUTINE TrilinosCleanup
121 !END INTERFACE
122 !#endif
123 
124  INTEGER :: runs
125 
126 
127  firstload = .true.
128  firsttime=.true.
129  initialize = 0
130 
131  ! Start the watches, store later
132  !--------------------------------
133  rt0 = realtime()
134  ct0 = cputime()
135 
136 
137  ! If parallel execution requested, initialize parallel environment:
138  !------------------------------------------------------------------
139  parallelenv => parallelinit()
140  outputpe = parenv % MyPE
141 
142 !tt = realtime()
143  IF ( firsttime ) THEN
144  !
145  ! Print banner to output:
146 #include "../config.h"
147  ! -----------------------
148  noargs = iargc()
149  ! Info Level is always true until the model has been read!
150  ! This makes it possible to cast something
151  silent = .false.
152  version = .false.
153  IF( noargs > 0 ) THEN
154  DO i = 1, noargs
155  CALL getarg( i,optionstring )
156  silent = silent .OR. &
157  ( optionstring=='-s' .OR. optionstring=='--silent' )
158  version = version .OR. &
159  ( optionstring=='-v' .OR. optionstring == '--version' )
160  END DO
161  END IF
162 
163  ! Set number of OpenMP threads
164  nthreads = 1
165  !$ nthreads = omp_get_max_threads()
166  IF (nthreads > 1) THEN
167  ! Check if OMP_NUM_THREADS environment variable is set
168  CALL envir( 'OMP_NUM_THREADS'//char(0), threads, tlen )
169  IF (tlen==0) THEN
170  CALL warn('MAIN','OMP_NUM_THREADS not set. Using only 1 thread.')
171  nthreads = 1
172  ! Set number of threads to 1
173  !$ CALL omp_set_num_threads(nthreads)
174 #ifdef HAVE_MKL
175  CALL mkl_set_num_threads(nthreads)
176 #endif
177  END IF
178  END IF
179 
180 
181  IF( .NOT. silent ) THEN
182  CALL info( 'MAIN', ' ')
183  CALL info( 'MAIN', '=============================================================')
184  CALL info( 'MAIN', 'MyElmerSolver finite element software, Welcome! ')
185  CALL info( 'MAIN', 'This program is free software licensed under (L)GPL ')
186  CALL info( 'MAIN', 'Copyright 1st April 1995 - , CSC - IT Center for Science Ltd.')
187  CALL info( 'MAIN', 'Webpage http://www.csc.fi/elmer, Email elmeradm@csc.fi ')
188  CALL info( 'MAIN', 'Library version: ' // version &
189 #ifdef REVISION
190  // ' (Rev: ' // revision // ')' )
191 #else
192  )
193 #endif
194  IF ( parenv % PEs > 1 ) &
195  CALL info( 'MAIN', ' Running in parallel using ' // &
196  trim(i2s(parenv % PEs)) // ' tasks.')
197 
198  ! Print out number of threads in use
199  IF ( nthreads > 1 ) &
200  CALL info('MAIN', ' Running in parallel with ' // &
201  trim(i2s(nthreads)) // ' threads per task.')
202 
203 #ifdef HAVE_HYPRE
204  CALL info( 'MAIN', ' HYPRE library linked in.')
205 #endif
206 #ifdef HAVE_TRILINOS
207  CALL info( 'MAIN', ' Trilinos library linked in.')
208 #endif
209 #ifdef HAVE_MUMPS
210  CALL info( 'MAIN', ' MUMPS library linked in.')
211 #endif
212 #ifdef HAVE_CHOLMOD
213  CALL info( 'MAIN', ' CHOLMOD library linked in.')
214 #endif
215 #ifdef HAVE_SUPERLU
216  CALL info( 'MAIN', ' SUPERLU library linked in.')
217 #endif
218 #ifdef HAVE_PARDISO
219  CALL info( 'MAIN', ' PARDISO library linked in.')
220 #endif
221 #ifdef HAVE_MKL
222  CALL info( 'MAIN', ' Intel MKL linked in.' )
223 #endif
224  CALL info( 'MAIN', '=============================================================')
225  END IF
226 
227  IF( version ) RETURN
228 
229  CALL initializeelementdescriptions
230 
231  WRITE(*,*) "Initialize: Setting FirstTime to FALSE"
232  firsttime = .false.
233  END IF
234 
235  ! Read input file name either as an argument, or from the default file:
236  !----------------------------------------------------------------------
237  gotmodelname = .false.
238  IF ( parenv % PEs <= 1 .AND. noargs > 0 ) THEN
239  CALL getarg( 1,modelname )
240  IF( modelname(1:1) /= '-') THEN
241  gotmodelname = .true.
242  IF ( noargs > 1 ) CALL getarg( 2,eq )
243  END IF
244  END IF
245 
246  IF( .NOT. gotmodelname ) THEN
247  OPEN( 1, file='ELMERSOLVER_STARTINFO', status='OLD', err=10 )
248  READ(1,'(a)') modelname
249  CLOSE(1)
250  END IF
251 
252 !------------------------------------------------------------------------------
253 ! Read element definition file, and initialize element types
254 !------------------------------------------------------------------------------
255  IF ( initialize==1 ) THEN
256  CALL freemodel(currentmodel)
257  firstload=.true.
258  END IF
259 
260 
261 !------------------------------------------------------------------------------
262 ! Read Model and mesh from Elmer mesh data base
263 !------------------------------------------------------------------------------
264 ! DO WHILE( .TRUE. )
265 
266  IF ( initialize==2 ) goto 1
267 
268  IF( myverbosity > 3) WRITE(*,*) 'FirstLoad = ', firstload
269 
270  IF ( firstload ) THEN
271  IF( .NOT. silent ) THEN
272  CALL info( 'MAIN', ' ')
273  CALL info( 'MAIN', ' ')
274  CALL info( 'MAIN', '-------------------------------------')
275  CALL info( 'MAIN', 'Reading Model: '//trim( modelname) )
276  END IF
277 
278  INQUIRE(unit=infileunit, opened=gotit)
279  IF ( gotit ) CLOSE(infileunit)
280 
281  IF( myverbosity > 3) WRITE(*,*) 'ModelName = ', modelname
282 
283  OPEN( unit=infileunit, action='Read',file=modelname,status='OLD',err=20 )
284  IF( myverbosity > 3) WRITE(*,*) 'Before LoadModel call'
285  currentmodel => loadmodel(modelname,.false.,parenv % PEs,parenv % MyPE )
286  IF( myverbosity > 3) WRITE(*,*) 'After LoadModel call'
287 
288  ! Optionally perform simple extrusion to increase the dimension of the mesh
289  !----------------------------------------------------------------------------------
290  extrudelevels=getinteger(currentmodel % Simulation,'Extruded Mesh Levels',found)
291  IF(extrudelevels>1) THEN
292  extrudedmesh => meshextrude(currentmodel % Meshes, extrudelevels-2)
293  DO i=1,currentmodel % NumberOfSolvers
294  IF(ASSOCIATED(currentmodel % Solvers(i) % Mesh,currentmodel % Meshes)) &
295  currentmodel % Solvers(i) % Mesh => extrudedmesh
296  END DO
297  extrudedmesh % Next => currentmodel % Meshes % Next
298  currentmodel % Meshes => extrudedmesh
299 
300  ! If periodic BC given, compute boundary mesh projector:
301  ! ------------------------------------------------------
302  DO i = 1,currentmodel % NumberOfBCs
303  IF(ASSOCIATED(currentmodel % Bcs(i) % PMatrix)) &
304  CALL freematrix( currentmodel % BCs(i) % PMatrix )
305  currentmodel % BCs(i) % PMatrix => null()
306  k = listgetinteger( currentmodel % BCs(i) % Values, 'Periodic BC', gotit )
307  IF( gotit ) THEN
308  currentmodel % BCs(i) % PMatrix => periodicprojector( currentmodel, extrudedmesh, i, k )
309  END IF
310  END DO
311  END IF
312 
313  ! If requested perform coordinate transformation directly after is has been obtained.
314  ! Don't maintain the original mesh.
315  !----------------------------------------------------------------------------------
316  coordtransform = listgetstring(currentmodel % Simulation,'Coordinate Transformation',gotit)
317  IF( gotit ) THEN
318  ! Masoud
319  WRITE(*,*) 'Initialize: Coordinate Transformation is activated for Elmer.\n'
320  ! Masoud : End
321  CALL coordinatetransformation( currentmodel % Meshes, coordtransform, &
322  currentmodel % Simulation, .true. )
323  END IF
324  ! IF( ListCheckPresent(CurrentModel % Simulation,'Coordinate Transformation') ) THEN
325  ! CALL CoordinateTransformation( CurrentModel % Meshes, .TRUE. )
326  ! END IF
327 
328  IF(.NOT. silent ) THEN
329  CALL info( 'MAIN', '-------------------------------------')
330  END IF
331  ELSE
332  IF ( initialize==3 ) THEN
333  INQUIRE(unit=infileunit, opened=gotit)
334  IF ( gotit ) CLOSE(infileunit)
335  OPEN( unit=infileunit, action='Read', &
336  file=modelname,status='OLD',err=20 )
337  END IF
338 
339  IF (myverbosity > 3) THEN
340  WRITE(6,*) 'Initialize: CurrentModel % NumberOfNodes&
341  = ', currentmodel % NumberOfNodes
342  END IF
343 
344 ! IF ( .NOT.ReloadInputFile(CurrentModel) ) EXIT
345 
346  mesh => currentmodel % Meshes
347  DO WHILE( ASSOCIATED(mesh) )
348  mesh % SavesDone = 0
349  mesh => mesh % Next
350  END DO
351  END IF
352 
353 1 CONTINUE
354 
355  CALL listaddlogical( currentmodel % Simulation, &
356  'Initialization Phase', .true. )
357 
358  ! Save the start time to the list so that it may be retrieved when necessary
359  ! This could perhaps also be a global variable etc, but this will do for now.
360  !-------------------------------------------------------------------------
361  IF( listgetlogical( currentmodel % Simulation,'Simulation Timing',gotit) ) THEN
362  CALL listaddconstreal( currentmodel % Simulation,'cputime0',ct0 )
363  CALL listaddconstreal( currentmodel % Simulation,'realtime0',rt0 )
364  END IF
365 #if 0
366 !------------------------------------------------------------------------------
367 ! Check for transient case
368 !------------------------------------------------------------------------------
369  eq = listgetstring( currentmodel % Simulation, 'Simulation Type', gotit )
370  scanning = eq == 'scanning'
371  transient = eq == 'transient'
372 
373 
374 !------------------------------------------------------------------------------
375 ! To more conveniently support the use of VTK based visualization there
376 ! is a hack that recognizes VTU suffix and creates a instance of output solver.
377 ! Note that this is really quite a dirty hack, and is not a good example.
378 !-----------------------------------------------------------------------------
380 
381 !------------------------------------------------------------------------------
382 ! Figure out what (flow,heat,stress,...) should be computed, and get
383 ! memory for the dofs
384 !------------------------------------------------------------------------------
385  CALL addsolvers()
386 
387 !------------------------------------------------------------------------------
388 ! Time integration and/or steady state steps
389 !------------------------------------------------------------------------------
390  IF ( transient .OR. scanning ) THEN
391  timesteps => listgetintegerarray( currentmodel % Simulation, &
392  'Timestep Intervals', gotit )
393 
394  IF ( .NOT.gotit ) THEN
395  CALL fatal( ' ', 'Keyword > Timestep Intervals < MUST be ' // &
396  'defined for transient and scanning simulations' )
397  END IF
398  timestepsizes => listgetconstrealarray( currentmodel % Simulation, &
399  'Timestep Sizes', gotit )
400 
401  IF ( .NOT.gotit ) THEN
402  IF( scanning .OR. listcheckpresent( currentmodel % Simulation,'Timestep Size') ) THEN
403  ALLOCATE(timestepsizes(SIZE(timesteps),1))
404  timestepsizes = 1.0_dp
405  ELSE
406  CALL fatal( ' ', 'Keyword [Timestep Sizes] MUST be ' // &
407  'defined for time dependent simulations' )
408  stop
409  END IF
410  END IF
411 
412  timeintervals = SIZE(timesteps)
413 
414  coupledmaxiter = listgetinteger( currentmodel % Simulation, &
415  'Steady State Max Iterations', gotit, minv=1 )
416  IF ( .NOT. gotit ) coupledmaxiter = 1
417 !------------------------------------------------------------------------------
418  ELSE
419 !------------------------------------------------------------------------------
420 ! Steady state
421 !------------------------------------------------------------------------------
422  ALLOCATE(timesteps(1))
423 
424  timesteps(1) = listgetinteger( currentmodel % Simulation, &
425  'Steady State Max Iterations', gotit,minv=1 )
426  IF ( .NOT. gotit ) timesteps(1)=1
427 
428  ALLOCATE(timestepsizes(1,1))
429  timestepsizes(1,1) = 1.0d0
430 
431  coupledmaxiter = 1
432  timeintervals = 1
433  END IF
434 
435  IF ( firstload ) &
436  ALLOCATE( stime(1), sstep(1), sinterval(1), ssize(1), &
437  steadyit(1), nonlinit(1), sprevsizes(1,5), speriodic(1) )
438 
439  dt = 0._dp
440 
441  stime = 0._dp
442  sstep = 0
443  speriodic = 0._dp
444 
445  ssize = dt
446  sprevsizes = 0_dp
447 
448  sinterval = 0._dp
449 
450  steadyit = 0
451  nonlinit = 0
452 
453  coupledminiter = listgetinteger( currentmodel % Simulation, &
454  'Steady State Min Iterations', gotit )
455 
456 !------------------------------------------------------------------------------
457 ! Add coordinates and simulation time to list of variables so that
458 ! coordinate dependent parameter computing routines can ask for
459 ! them...
460 !------------------------------------------------------------------------------
461  IF ( firstload ) CALL addmeshcoordinatesandtime
462 
463 !------------------------------------------------------------------------------
464 ! Get Output File Options
465 !------------------------------------------------------------------------------
466 
467  outputintervals => listgetintegerarray( currentmodel % Simulation, &
468  'Output Intervals', gotit )
469  IF ( .NOT. gotit ) THEN
470  ALLOCATE( outputintervals(SIZE(timesteps)) )
471  outputintervals = 1
472  END IF
473 
474 
475  ! Initial Conditions:
476  ! -------------------
477  IF ( firstload ) CALL setinitialconditions
478 
479 
480  ! Compute the total number of steps that will be saved to the files
481  ! Particularly look if the last step will be saved, or if it has
482  ! to be saved separately.
483  !------------------------------------------------------------------
484  totaltimesteps = 0
485  lastsaved = .true.
486  DO interval=1,timeintervals
487  DO timestep = 1,timesteps(interval)
488  IF( outputintervals(interval) == 0 ) cycle
489  lastsaved = .false.
490  IF ( mod(timestep-1, outputintervals(interval))==0 ) THEN
491  lastsaved = .true.
492  totaltimesteps = totaltimesteps + 1
493  END IF
494  END DO
495  END DO
496 
497  DO i=1,currentmodel % NumberOfSolvers
498  solver => currentmodel % Solvers(i)
499  IF(.NOT.ASSOCIATED(solver % Variable)) cycle
500  IF(.NOT.ASSOCIATED(solver % Variable % Values)) cycle
501  when = listgetstring( solver % Values, 'Exec Solver', gotit )
502  IF ( gotit ) THEN
503  IF ( when == 'after simulation' .OR. when == 'after all' ) THEN
504  lastsaved = .false.
505  END IF
506  ELSE
507  IF ( solver % SolverExecWhen == solver_exec_after_all ) THEN
508  lastsaved = .false.
509  END IF
510  END IF
511  END DO
512 
513  IF ( .NOT.lastsaved ) totaltimesteps = totaltimesteps + 1
514  IF( totaltimesteps == 0 ) totaltimesteps = 1
515 
516  CALL listaddlogical( currentmodel % Simulation, &
517  'Initialization Phase', .false. )
518 
519  firstload = .false.
520  IF (myverbosity > 3) THEN
521  WRITE(6,*) 'MyElmerSolver 493 Initialize = ', initialize
522  END IF
523 ! IF ( Initialize == 1 ) EXIT
524 #endif
525 
526 ! END DO
527 
528  IF (myverbosity > 3) THEN
529  WRITE(6,*) 'MyRun: CurrentModel % NumberOfNodes = ', &
530  currentmodel % NumberOfNodes
531  END IF
532 
533 !------------------------------------------------------------------------------
534 ! Check for transient case
535 !------------------------------------------------------------------------------
536  eq = listgetstring( currentmodel % Simulation, 'Simulation Type', gotit )
537  scanning = eq == 'scanning'
538  transient = eq == 'transient'
539 
540 
541 !------------------------------------------------------------------------------
542 ! To more conveniently support the use of VTK based visualization there
543 ! is a hack that recognizes VTU suffix and creates a instance of output solver.
544 ! Note that this is really quite a dirty hack, and is not a good example.
545 !-----------------------------------------------------------------------------
546  IF (myverbosity > 3) THEN
547  write(6,*) 'calling AddVtuOutputSolverHack'
548  END IF
550  IF (myverbosity > 3) THEN
551  write(6,*) 'leaving AddVtuOutputSolverHack'
552  END IF
553 
554 !------------------------------------------------------------------------------
555 ! Figure out what (flow,heat,stress,...) should be computed, and get
556 ! memory for the dofs
557 !------------------------------------------------------------------------------
558  IF (myverbosity > 3) THEN
559  write(6,*) 'calling AddSolvers'
560  END IF
561  CALL addsolvers()
562  IF (myverbosity > 3) THEN
563  write(6,*) 'leaving AddSolvers'
564  END IF
565 
566 !------------------------------------------------------------------------------
567 ! Time integration and/or steady state steps
568 !------------------------------------------------------------------------------
569  IF ( transient .OR. scanning ) THEN
570  timesteps => listgetintegerarray( currentmodel % Simulation, &
571  'Timestep Intervals', gotit )
572 
573  IF ( .NOT.gotit ) THEN
574  CALL fatal( ' ', 'Keyword > Timestep Intervals < MUST be ' // &
575  'defined for transient and scanning simulations' )
576  END IF
577  timestepsizes => listgetconstrealarray( currentmodel % Simulation, &
578  'Timestep Sizes', gotit )
579 
580  IF ( .NOT.gotit ) THEN
581  IF( scanning .OR. listcheckpresent( currentmodel % Simulation,'Timestep Size') ) THEN
582  ALLOCATE(timestepsizes(SIZE(timesteps),1))
583  timestepsizes = 1.0_dp
584  ELSE
585  CALL fatal( ' ', 'Keyword [Timestep Sizes] MUST be ' // &
586  'defined for time dependent simulations' )
587  stop
588  END IF
589  END IF
590 
591  timeintervals = SIZE(timesteps)
592 
593  coupledmaxiter = listgetinteger( currentmodel % Simulation, &
594  'Steady State Max Iterations', gotit, minv=1 )
595  IF ( .NOT. gotit ) coupledmaxiter = 1
596 !------------------------------------------------------------------------------
597  ELSE
598 !------------------------------------------------------------------------------
599 ! Steady state
600 !------------------------------------------------------------------------------
601  ALLOCATE(timesteps(1))
602 
603  timesteps(1) = listgetinteger( currentmodel % Simulation, &
604  'Steady State Max Iterations', gotit,minv=1 )
605  IF ( .NOT. gotit ) timesteps(1)=1
606 
607  ALLOCATE(timestepsizes(1,1))
608  timestepsizes(1,1) = 1.0d0
609 
610  coupledmaxiter = 1
611  timeintervals = 1
612  END IF
613 
614  IF ( firstload ) &
615  ALLOCATE( stime(1), sstep(1), sinterval(1), ssize(1), &
616  steadyit(1), nonlinit(1), sprevsizes(1,5), speriodic(1) )
617 
618  dt = 0._dp
619 
620  stime = 0._dp
621  sstep = 0
622  speriodic = 0._dp
623 
624  ssize = dt
625  sprevsizes = 0_dp
626 
627  sinterval = 0._dp
628 
629  steadyit = 0
630  nonlinit = 0
631 
632  coupledminiter = listgetinteger( currentmodel % Simulation, &
633  'Steady State Min Iterations', gotit )
634 
635 !------------------------------------------------------------------------------
636 ! Add coordinates and simulation time to list of variables so that
637 ! coordinate dependent parameter computing routines can ask for
638 ! them...
639 !------------------------------------------------------------------------------
640  IF ( firstload ) CALL addmeshcoordinatesandtime
641 
642 !------------------------------------------------------------------------------
643 ! Get Output File Options
644 !------------------------------------------------------------------------------
645 
646  outputintervals => listgetintegerarray( currentmodel % Simulation, &
647  'Output Intervals', gotit )
648 
649  IF ( .NOT. gotit ) THEN
650  ALLOCATE( outputintervals(SIZE(timesteps)) )
651  outputintervals = 1
652  END IF
653 
654  IF (myverbosity > 3) THEN
655  write(6,*) 'MyRunI SIZE(OutputIntervals) =', SIZE(outputintervals)
656  write(6,*) 'MyRunI OutputIntervals(1) = ',outputintervals(1)
657  END IF
658 
659  ! Initial Conditions:
660  ! -------------------
661  IF ( firstload ) CALL setinitialconditions
662 
663 
664  ! Compute the total number of steps that will be saved to the files
665  ! Particularly look if the last step will be saved, or if it has
666  ! to be saved separately.
667  !------------------------------------------------------------------
668  totaltimesteps = 0
669  lastsaved = .true.
670  DO interval=1,timeintervals
671  DO timestep = 1,timesteps(interval)
672  IF( outputintervals(interval) == 0 ) cycle
673  lastsaved = .false.
674  IF ( mod(timestep-1, outputintervals(interval))==0 ) THEN
675  lastsaved = .true.
676  totaltimesteps = totaltimesteps + 1
677  END IF
678  END DO
679  END DO
680 
681  DO i=1,currentmodel % NumberOfSolvers
682  solver => currentmodel % Solvers(i)
683  IF(.NOT.ASSOCIATED(solver % Variable)) cycle
684  IF(.NOT.ASSOCIATED(solver % Variable % Values)) cycle
685  when = listgetstring( solver % Values, 'Exec Solver', gotit )
686  IF ( gotit ) THEN
687  IF ( when == 'after simulation' .OR. when == 'after all' ) THEN
688  lastsaved = .false.
689  END IF
690  ELSE
691  IF ( solver % SolverExecWhen == solver_exec_after_all ) THEN
692  lastsaved = .false.
693  END IF
694  END IF
695  END DO
696 
697  IF ( .NOT.lastsaved ) totaltimesteps = totaltimesteps + 1
698  IF( totaltimesteps == 0 ) totaltimesteps = 1
699 
700  CALL listaddlogical( currentmodel % Simulation, &
701  'Initialization Phase', .false. )
702 
703  firstload = .false.
704  IF (myverbosity > 3) THEN
705  WRITE(6,*) 'MyElmerSolver 493 Initialize = ', initialize
706  write(6,*) 'TimeIntervals from inside RunI = ', timeintervals
707  WRITE(*,*) 'NumberOfSolvers = ', currentmodel % NumberOfSolvers
708  END IF
709 
710 ! IF ( Initialize == 1 ) EXIT
711 
712  DO i=1,currentmodel % NumberOfSolvers
713  solver => currentmodel % Solvers(i)
714  IF (myverbosity > 3) THEN
715  WRITE(*,*) 'Solver % PROCEDURE = ', solver % PROCEDURE
716  END IF
717  IF ( solver % PROCEDURE==0 ) cycle
718  IF ( solver % SolverExecWhen == solver_exec_ahead_all ) THEN
719  IF (myverbosity > 3) THEN
720  WRITE(*,*) 'Calling SolverActivate for Solver', i
721  END IF
722  CALL solveractivate( currentmodel,solver,dt,transient )
723  END IF
724  END DO
725 
726  runs = runs + 1
727 
728  RETURN
729 
730 10 CONTINUE
731  CALL fatal( 'MyElmerSolver', 'Unable to find ELMERSOLVER_STARTINFO, can not execute.' )
732 
733 20 CONTINUE
734  CALL fatal( 'MyElmerSolver', 'Unable to find input file [' // &
735  trim(modelname) // '], can not execute.' )
736 
737 
738 !------------------------------------------------------------------------------
739  END SUBROUTINE elmerinitialize
740 !------------------------------------------------------------------------------
741 
subroutine setinitialconditions()
Sets initial conditions for the fields.
subroutine addvtuoutputsolverhack()
subroutine addsolvers()
Adds flags for active solvers.
subroutine addmeshcoordinatesandtime()
Adds coordinate and time variables to the current mesh structure.
subroutine elmerinitialize(runs)
ElmerLib Elmer library routines
Definition: Initialize.F90:55
ElmerLib Elmer library routines