ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
initializeFlow.F90
Go to the documentation of this file.
2 
3  use constants, only: inttype, realtype, maxstringlen
4 
5  implicit none
6  save
7 
8 contains
9 
10  subroutine referencestate
11  !
12  ! The original version has been nuked since the computations are
13  ! no longer necessary when calling from python
14  ! This is the most compliclated routine in all of ADflow. It is
15  ! stupidly complicated. This is most likely the reason your
16  ! derivatives are wrong. You don't understand this routine
17  ! and its effects.
18  ! This routine *requries* the following as input:
19  ! Mach, pInfDim, TInfDim, rhoInfDim, rGasDim (machCoef non-SA
20  ! turbulence only)
21  ! Optionally, pRef, rhoRef and Tref are used if they are
22  ! are non-negative. This only happens when you want the equations
23  ! normalized by values other than the freestream
24  ! * This routine computes as output:
25  ! * muInfDim, (unused anywhere in code)
26  ! pRef, rhoRef, Tref, muRef, timeRef ('dimensional' reference)
27  ! pInf, pInfCorr, rhoInf, uInf, rGas, muInf, gammaInf and wInf
28  ! (Non-dimensionalized values used in actual computations)
29  !
30  use constants
31  use paramturb
32  use inputphysics, only: equations, mach, machcoef, &
36  muinfdim, &
39  nw, nwf, kpresent, winf
40  use flowutils, only: computegamma, etot
41  use turbutils, only: sanuknowneddyratio
42  implicit none
43 
44  integer(kind=intType) :: sps, nn, mm, ierr
45  real(kind=realtype) :: gm1, ratio
46  real(kind=realtype) :: nuinf, ktmp, uinf2
47  real(kind=realtype) :: vinf, zinf, tmp1(1), tmp2(1)
48 
49 
50  ! Compute the dimensional viscosity from Sutherland's law
52  * ((tsuthdim + ssuthdim) / (tinfdim + ssuthdim)) &
53  * ((tinfdim / tsuthdim)**1.5_realtype)
54 
55  ! Set the reference values. They *COULD* be different from the
56  ! free-stream values for an internal flow simulation. For now,
57  ! we just use the actual free stream values.
58  pref = pinfdim
59  tref = tinfdim
61 
62  ! Compute the value of muRef, such that the nonDimensional
63  ! equations are identical to the dimensional ones.
64  ! Note that in the non-dimensionalization of muRef there is
65  ! a reference length. However this reference length is 1.0
66  ! in this code, because the coordinates are converted to
67  ! meters.
68 
69  muref = sqrt(pref * rhoref)
70 
71  ! Compute timeRef for a correct nonDimensionalization of the
72  ! unsteady equations. Some story as for the reference viscosity
73  ! concerning the reference length.
74 
75  timeref = sqrt(rhoref / pref)
76  href = pref / rhoref
77  uref = sqrt(href)
78 
79  ! Compute the nonDimensional pressure, density, velocity,
80  ! viscosity and gas constant.
81 
82  pinf = pinfdim / pref
84  uinf = mach * sqrt(gammainf * pinf / rhoinf)
85  rgas = rgasdim * rhoref * tref / pref
86  muinf = muinfdim / muref
87  tmp1(1) = tinfdim
88  call computegamma(tmp1, tmp2, 1)
89  gammainf = tmp2(1)
90 
91  ! ----------------------------------------
92  ! Compute the final wInf
93  ! ----------------------------------------
94 
95  ! Allocate the memory for wInf if necessary
96 #ifndef USE_TAPENADE
97  if (allocated(winf)) deallocate (winf)
98  allocate (winf(nw), stat=ierr)
99 #endif
100 
101  ! zero out the winf first
102  winf(:) = zero
103 
104  ! Set the reference value of the flow variables, except the total
105  ! energy. This will be computed at the end of this routine.
106 
107  winf(irho) = rhoinf
108  winf(ivx) = uinf * veldirfreestream(1)
109  winf(ivy) = uinf * veldirfreestream(2)
110  winf(ivz) = uinf * veldirfreestream(3)
111 
112  ! Compute the velocity squared based on MachCoef. This gives a
113  ! better indication of the 'speed' of the flow so the turubulence
114  ! intensity ration is more meaningful especially for moving
115  ! geometries. (Not used in SA model)
116 
117  uinf2 = machcoef * machcoef * gammainf * pinf / rhoinf
118 
119  ! Set the turbulent variables if transport variables are to be
120  ! solved. We should be checking for RANS equations here,
121  ! however, this code is included in block res. The issue is
122  ! that for frozen turbulence (or ANK jacobian) we call the
123  ! block_res with equationType set to Laminar even though we are
124  ! actually solving the rans equations. The issue is that, the
125  ! freestream turb variables will be changed to zero, thus
126  ! changing the solution. Insteady we check if nw > nwf which
127  ! will accomplish the same thing.
128 
129  if (nw > nwf) then
130 
131  nuinf = muinf / rhoinf
132 
133  select case (turbmodel)
134 
136 
138 
139  !=============================================================
140 
142 
143  winf(itu1) = 1.5_realtype * uinf2 * turbintensityinf**2
144  winf(itu2) = winf(itu1) / (eddyvisinfratio * nuinf)
145 
146  !=============================================================
147 
148  case (ktau)
149 
150  winf(itu1) = 1.5_realtype * uinf2 * turbintensityinf**2
151  winf(itu2) = eddyvisinfratio * nuinf / winf(itu1)
152 
153  !=============================================================
154 
155  case (v2f)
156 
157  winf(itu1) = 1.5_realtype * uinf2 * turbintensityinf**2
158  winf(itu2) = 0.09_realtype * winf(itu1)**2 &
159  / (eddyvisinfratio * nuinf)
160  winf(itu3) = 0.666666_realtype * winf(itu1)
161  winf(itu4) = 0.0_realtype
162 
163  end select
164 
165  end if
166 
167  ! Set the value of pInfCorr. In case a k-equation is present
168  ! add 2/3 times rho*k.
169 
170  pinfcorr = pinf
171  if (kpresent) pinfcorr = pinf + two * third * rhoinf * winf(itu1)
172 
173  ! Compute the free stream total energy.
174 
175  ktmp = zero
176  if (kpresent) ktmp = winf(itu1)
177  vinf = zero
178  zinf = zero
179  call etot(rhoinf, uinf, vinf, zinf, pinfcorr, ktmp, &
180  winf(irhoe), kpresent)
181 
182  end subroutine referencestate
183 
184  ! ----------------------------------------------------------------------
185  ! |
186  ! No Tapenade Routine below this line |
187  ! |
188  ! ----------------------------------------------------------------------
189 
190 #ifndef USE_TAPENADE
191  subroutine infchangecorrection(oldWinf, correctionTol, correctionType)
192  ! Adjust the flow states to a change in wInf
193  use constants
194  use blockpointers, only: il, jl, kl, w, ndom, d2wall
195  use flowvarrefstate, only: winf, nwf, nw
196  use inputphysics, only: equations
198  use haloexchange, only: whalo2
199  use flowutils, only: adjustinflowangle
200  use oversetdata, only: oversetpresent
201  use iteration, only: currentlevel
203  use bcroutines, only: applyallbc_block
205  use communication, only: myid
206  implicit none
207 
208  real(kind=realtype), intent(in), dimension(nwf) :: oldwinf
209  real(kind=realtype), intent(in) :: correctiontol
210  character*(*), intent(in) :: correctionType
211  integer(kind=intType) :: sps, nn, i, j, k, l
212  real(kind=realtype) :: deltawinf(nwf)
213 
214  ! variables for rotation update
215  real(kind=realtype), dimension(3) :: vec1, vec2, vcell
216  real(kind=realtype), dimension(3, 3) :: rotmat
217  real(kind=realtype) :: mag1, mag2
218 
219  ! Make sure we have the updated wInf
220  call adjustinflowangle()
221  call referencestate
222 
223  ! first check if the state changed enough to warrant a correction. This tolerance
224  ! can be adjusted from python. Default is 1e-12.
225  deltawinf = winf(1:nwf) - oldwinf(1:nwf)
226 
227  if (mynorm2(deltawinf) < correctiontol) then
228  ! The change deltaWinf is so small, (or zero) don't do the
229  ! update and just return. This will save some time when the
230  ! solver is called with the same AP conditions multiple times,
231  ! such as during a GS AS solution
232  return
233  end if
234 
235  ! The state changed enough so we will run a correction update.
236  ! Check which correction type we are using
237  if (correctiontype .eq. "offset") then
238 
239  ! Loop over all the blocks, adding the offset between the old and new Winf
240  do sps = 1, ntimeintervalsspectral
241  do nn = 1, ndom
242  call setpointers(nn, 1_inttype, sps)
243  do k = 2, kl
244  do j = 2, jl
245  do i = 2, il
246  do l = 1, nwf
247  w(i, j, k, l) = w(i, j, k, l) + deltawinf(l)
248  end do
249  end do
250  end do
251  end do
252  call applyallbc_block(.true.)
253  end do
254  end do
255 
256  else if (correctiontype .eq. "rotate") then
257 
258  ! Rotate the flow velocities by the rotation change in Winf.
259  ! We also scale the velocity magnitudes based on the magnitude change.
260  ! vec1 is the old free stream velocity vector, vec2 is the new one.
261  vec1 = oldwinf(ivx:ivz)
262  vec2 = winf(ivx:ivz)
263 
264  ! get the general rotation matrix
265  call getrotmatrix(vec1, vec2, rotmat)
266 
267  ! Compute the magnitude of vec1 and vec2 so that we can adjust the magnitude of the velocity
268  mag1 = mynorm2(vec1)
269  mag2 = mynorm2(vec2)
270 
271  ! Loop over all the blocks, rotate each cell's velocity
272  do sps = 1, ntimeintervalsspectral
273  do nn = 1, ndom
274  call setpointers(nn, 1_inttype, sps)
275  do k = 2, kl
276  do j = 2, jl
277  do i = 2, il
278  ! pick out the cell velocity
279  vcell = w(i, j, k, ivx:ivz)
280  ! rotate and put back
281  w(i, j, k, ivx:ivz) = matmul(rotmat, vcell)
282  ! scale the magnitudes
283  w(i, j, k, ivx:ivz) = w(i, j, k, ivx:ivz) * mag2 / mag1
284  ! we also add the offset to the density and energy
285  w(i, j, k, irho) = w(i, j, k, irho) + deltawinf(irho)
286  w(i, j, k, irhoe) = w(i, j, k, irhoe) + deltawinf(irhoe)
287  end do
288  end do
289  end do
290  call applyallbc_block(.true.)
291  end do
292  end do
293 
294  end if
295 
296  ! Exchange values
297  call whalo2(currentlevel, 1_inttype, nw, .true., .true., .true.)
298 
299  ! Need to re-apply the BCs. The reason is that BC halos behind
300  ! interpolated cells need to be recomputed with their new
301  ! interpolated values from actual compute cells. Only needed for
302  ! overset.
303  if (oversetpresent) then
304  do sps = 1, ntimeintervalsspectral
305  do nn = 1, ndom
306  call setpointers(nn, 1, sps)
307  if (equations == ransequations) then
308  call bcturbtreatment
309  call applyallturbbcthisblock(.true.)
310  end if
311  call applyallbc_block(.true.)
312  end do
313  end do
314  end if
315 
316  end subroutine infchangecorrection
317 
318  ! Section out the BCdata setup so that it can by called from python when needed
320  ! sets the prescribed boundary data from the CGNS arrays
321 
322  use constants
323  use iteration, only: groundlevel
325  implicit none
326 
327  ! Allocate the memory for the prescribed boundary data at the
328  ! boundary faces and determine the data for both the fine grid.
329 
330  groundlevel = 1
331 
332  ! Determine the reference state.
333  call referencestate
334 
335  call setbcdatafinegrid(.true.)
336 
337  ! Determine the prescribed data on the coarse grid levels
338  ! by interpolation.
339 #ifndef USE_TAPENADE
341 #endif
342 
343  end subroutine updatebcdataalllevels
344 
345  subroutine initflow
346  !
347  ! initFlow allocates the
348  ! memory for and initializes the flow variables. In case a
349  ! restart is performed the owned variables are read from the
350  ! previous solution file(s).
351  !
352  use constants
353  use block, only: flowdoms
355  use variablereading, only: halosread
356 
357  implicit none
358  !
359  ! Local variables.
360  !
361  integer :: ierr
362 
363  integer(kind=intType) :: sps, level, nLevels
364 
365  ! Determine the number of multigrid levels.
366 
367  nlevels = ubound(flowdoms, 2)
368 
369  ! As some boundary conditions can be treated in multiple ways,
370  ! some memory allocated must be released again.
371 
372  call releaseextramembcs
373 
374  ! Determine for the time spectral mode the matrices for the
375  ! time derivatives.
377 
378  ! Loop over the number of spectral solutions to allocate
379  ! the memory for the w-variables and p on the fine grid.
380  do sps = 1, ntimeintervalsspectral
381  call allocmemflovarpart1(sps, 1_inttype)
382  end do
383 
384  ! Allocate the memory for the solution variables on the coarse
385  ! grid levels and the memory for the dependent flow variables,
386  ! residuals, etc, on all multigrid levels.
387  do sps = 1, ntimeintervalsspectral
388  call allocmemflovarpart2(sps, 1_inttype)
389 
390  do level = 2, nlevels
391  call allocmemflovarpart1(sps, level)
392  call allocmemflovarpart2(sps, level)
393  end do
394  end do
395 
396  ! Initialize free stream field
397  call initflowfield
398 
399  ! Initialize the dependent flow variables and the halo values.
400 
402 
403  end subroutine initflow
404 
405  subroutine allocmemflovarpart1(sps, level)
406  !
407  ! allocMemFlovarPart1 allocates the memory for the flow
408  ! variables w and p for all the blocks on the given multigrid
409  ! level and spectral solution sps.
410  !
411  use constants
412  use block, only: flowdoms, ndom
413  use flowvarrefstate, only: nw, nwf, nt1, nt2
417  use iteration, only: noldlevels
418  use utils, only: terminate
419  implicit none
420  !
421  ! Subroutine arguments.
422  !
423  integer(kind=intType), intent(in) :: sps, level
424  !
425  ! Local variables.
426  !
427  integer :: ierr
428 
429  integer(kind=intType) :: nn
430  integer(kind=intType) :: il, jl, kl, ie, je, ke, ib, jb, kb
431 
432  ! Loop over the domains.
433 
434  domains: do nn = 1, ndom
435 
436  ! Store some dimensions a bit easier.
437 
438  il = flowdoms(nn, level, sps)%il
439  jl = flowdoms(nn, level, sps)%jl
440  kl = flowdoms(nn, level, sps)%kl
441 
442  ie = flowdoms(nn, level, sps)%ie
443  je = flowdoms(nn, level, sps)%je
444  ke = flowdoms(nn, level, sps)%ke
445 
446  ib = flowdoms(nn, level, sps)%ib
447  jb = flowdoms(nn, level, sps)%jb
448  kb = flowdoms(nn, level, sps)%kb
449 
450  ! Allocate the memory for the independent variables.
451  ! Memory is allocated for the turbulent variables (if any) if
452  ! the current level is smaller or equal to the multigrid start
453  ! level or if the turbulent transport equations are solved in
454  ! a coupled manner.
455 
456  if (level <= mgstartlevel .or. turbtreatment == coupled) then
457  allocate (flowdoms(nn, level, sps)%w(0:ib, 0:jb, 0:kb, 1:nw), &
458  stat=ierr)
459  else
460  allocate (flowdoms(nn, level, sps)%w(0:ib, 0:jb, 0:kb, 1:nwf), &
461  stat=ierr)
462  end if
463  if (ierr /= 0) &
464  call terminate("allocMemFlovarPart1", &
465  "Memory allocation failure for w")
466 
467  ! Alloc mem for nodal gradients
468  allocate (flowdoms(nn, level, sps)%ux(il, jl, kl), stat=ierr)
469  allocate (flowdoms(nn, level, sps)%uy(il, jl, kl), stat=ierr)
470  allocate (flowdoms(nn, level, sps)%uz(il, jl, kl), stat=ierr)
471 
472  allocate (flowdoms(nn, level, sps)%vx(il, jl, kl), stat=ierr)
473  allocate (flowdoms(nn, level, sps)%vy(il, jl, kl), stat=ierr)
474  allocate (flowdoms(nn, level, sps)%vz(il, jl, kl), stat=ierr)
475 
476  allocate (flowdoms(nn, level, sps)%wx(il, jl, kl), stat=ierr)
477  allocate (flowdoms(nn, level, sps)%wy(il, jl, kl), stat=ierr)
478  allocate (flowdoms(nn, level, sps)%wz(il, jl, kl), stat=ierr)
479 
480  allocate (flowdoms(nn, level, sps)%qx(il, jl, kl), stat=ierr)
481  allocate (flowdoms(nn, level, sps)%qy(il, jl, kl), stat=ierr)
482  allocate (flowdoms(nn, level, sps)%qz(il, jl, kl), stat=ierr)
483 
484  ! Allocate memory for the pressure.
485  allocate (flowdoms(nn, level, sps)%p(0:ib, 0:jb, 0:kb), stat=ierr)
486  if (ierr /= 0) &
487  call terminate("allocMemFlovarPart1", &
488  "Memory allocation failure for p")
489 
490  ! Allocate memory for the speed of sound squared
491  allocate (flowdoms(nn, level, sps)%aa(0:ib, 0:jb, 0:kb), stat=ierr)
492  if (ierr /= 0) &
493  call terminate("allocMemFlovarPart1", &
494  "Memory allocation failure for p")
495 
496  ! The eddy viscosity for eddy viscosity models.
497  ! Although a dependent variable, it is allocated on all grid
498  ! levels, because the eddy viscosity might be frozen in the
499  ! multigrid.
500 
501  ! Always allocate rev due to reverse mode AD- Peter Lyu. Also
502  ! zero so that it doesn't affect laminar cases.
503  allocate (flowdoms(nn, level, sps)%rev(0:ib, 0:jb, 0:kb), &
504  stat=ierr)
505  flowdoms(nn, level, sps)%rev = zero
506  if (ierr /= 0) &
507  call terminate("allocMemFlovarPart1", &
508  "Memory allocation failure for rev")
509  !endif
510 
511  ! If this is the finest grid some more memory must be allocated.
512 
513  fineleveltest: if (level == 1) then
514 
515  ! Allocate the memory for gamma and initialize it to
516  ! the constant gamma value.
517 
518  allocate (flowdoms(nn, level, sps)%gamma(0:ib, 0:jb, 0:kb), &
519  stat=ierr)
520  if (ierr /= 0) &
521  call terminate("allocMemFlovarPart1", &
522  "Memory allocation failure for gamma.")
523 
524  flowdoms(nn, level, sps)%gamma = gammaconstant
525 
526  ! The laminar viscosity for viscous computations.
527  ! Always allocate rlv due to reverse mode - Peter Lyu
528  !if( viscous ) then
529  allocate (flowdoms(nn, level, sps)%rlv(0:ib, 0:jb, 0:kb), &
530  stat=ierr)
531  if (ierr /= 0) &
532  call terminate("allocMemFlovarPart1", &
533  "Memory allocation failure for rlv")
534  !endif
535 
536  ! The state vectors in the past for unsteady computations.
537 
538  if (equationmode == unsteady .and. &
539  timeintegrationscheme == bdf) then
540  allocate ( &
541  flowdoms(nn, level, sps)%wOld(noldlevels, 2:il, 2:jl, 2:kl, nw), &
542  stat=ierr)
543  if (ierr /= 0) &
544  call terminate("allocMemFlovarPart1", &
545  "Memory allocation failure for wOld")
546 
547  ! Initialize wOld to zero, such that it is initialized.
548  ! The actual values do not matter.
549 
550  flowdoms(nn, level, sps)%wOld = zero
551 
552  ! Added by HDN
553  else if (equationmode == unsteady .and. &
554  timeintegrationscheme == md) then
555  allocate ( &
556  flowdoms(nn, level, sps)%wOld(noldlevels, 2:il, 2:jl, 2:kl, nw), &
557  stat=ierr)
558  if (ierr /= 0) &
559  call terminate("allocMemFlovarPart1", &
560  "Memory allocation failure for wOld")
561 
562  flowdoms(nn, level, sps)%wOld = zero
563  end if
564 
565  ! If this is the 1st spectral solution (note that we are
566  ! already on the finest grid) and the rans equations are
567  ! solved, allocate the memory for the arrays used for the
568  ! implicit boundary condition treatment. Normally this should
569  ! only be allocated for RANS but the derivative calcs require
570  ! these be allocated.
571 
572  sps1ranstest: if (sps == 1) then
573  allocate (flowdoms(nn, level, sps)%bmti1(je, ke, nt1:nt2, nt1:nt2), &
574  flowdoms(nn, level, sps)%bmti2(je, ke, nt1:nt2, nt1:nt2), &
575  flowdoms(nn, level, sps)%bmtj1(ie, ke, nt1:nt2, nt1:nt2), &
576  flowdoms(nn, level, sps)%bmtj2(ie, ke, nt1:nt2, nt1:nt2), &
577  flowdoms(nn, level, sps)%bmtk1(ie, je, nt1:nt2, nt1:nt2), &
578  flowdoms(nn, level, sps)%bmtk2(ie, je, nt1:nt2, nt1:nt2), &
579  flowdoms(nn, level, sps)%bvti1(je, ke, nt1:nt2), &
580  flowdoms(nn, level, sps)%bvti2(je, ke, nt1:nt2), &
581  flowdoms(nn, level, sps)%bvtj1(ie, ke, nt1:nt2), &
582  flowdoms(nn, level, sps)%bvtj2(ie, ke, nt1:nt2), &
583  flowdoms(nn, level, sps)%bvtk1(ie, je, nt1:nt2), &
584  flowdoms(nn, level, sps)%bvtk2(ie, je, nt1:nt2), &
585  stat=ierr)
586  if (ierr /= 0) &
587  call terminate("allocMemFlovarPart1", &
588  "Memory allocation failure for bmti1, etc")
589 
590  end if sps1ranstest
591  end if fineleveltest
592 
593  end do domains
594 
595  end subroutine allocmemflovarpart1
596 
597  ! ==================================================================
598 
599  subroutine allocmemflovarpart2(sps, level)
600  !
601  ! AllocMemFlovarPart2 allocates the memory for the dependent
602  ! flow variables and iteration variables for all the blocks on
603  ! the given multigrid level and spectral solution sps. Some
604  ! variables are only allocated on the coarser grids, e.g. the
605  ! multigrid forcing terms and the state vector upon entrance on
606  ! the mg level. Other variables are only allocated on the finest
607  ! mesh. These are typically dependent variables like laminar
608  ! viscosity, or residuals, time step, etc. Exceptions are
609  ! pressure and eddy viscosity. Although these are dependent
610  ! variables, they are allocated on all grid levels.
611  !
612  use block
613  use constants
614  use flowvarrefstate
615  use inputphysics
617  use inputiteration
618  use inputunsteady
619  use iteration
620  use utils, only: terminate
621  implicit none
622  !
623  ! Subroutine arguments.
624  !
625  integer(kind=intType), intent(in) :: sps, level
626  !
627  ! Local variables.
628  !
629  integer :: ierr
630 
631  integer(kind=intType) :: nn, mm
632  integer(kind=intType) :: il, jl, kl, ie, je, ke, ib, jb, kb
633 
634  ! Loop over the domains.
635 
636  domains: do nn = 1, ndom
637 
638  ! Store some dimensions a bit easier.
639 
640  il = flowdoms(nn, level, sps)%il
641  jl = flowdoms(nn, level, sps)%jl
642  kl = flowdoms(nn, level, sps)%kl
643 
644  ie = flowdoms(nn, level, sps)%ie
645  je = flowdoms(nn, level, sps)%je
646  ke = flowdoms(nn, level, sps)%ke
647 
648  ib = flowdoms(nn, level, sps)%ib
649  jb = flowdoms(nn, level, sps)%jb
650  kb = flowdoms(nn, level, sps)%kb
651 
652  ! Block is moving. Allocate the memory for s, sFaceI,
653  ! sFaceJ and sFaceK.
654 
655  allocate (flowdoms(nn, level, sps)%s(ie, je, ke, 3), &
656  flowdoms(nn, level, sps)%sFaceI(0:ie, je, ke), &
657  flowdoms(nn, level, sps)%sFaceJ(ie, 0:je, ke), &
658  flowdoms(nn, level, sps)%sFaceK(ie, je, 0:ke), stat=ierr)
659  if (ierr /= 0) &
660  call terminate("allocMemFlovarPart2", &
661  "Memory allocation failure for s, &
662  &sFaceI, sFaceJ and sFaceK.")
663 
664  ! Extra face velocities for ALE
665  if (equationmode == unsteady .and. useale) then
666  allocate ( &
667  flowdoms(nn, level, sps)%sVeloIALE(0:ie, je, ke, 3), &
668  flowdoms(nn, level, sps)%sVeloJALE(ie, 0:je, ke, 3), &
669  flowdoms(nn, level, sps)%sVeloKALE(ie, je, 0:ke, 3), &
670  flowdoms(nn, level, sps)%sFaceIALE(0:nalesteps, 0:ie, je, ke), &
671  flowdoms(nn, level, sps)%sFaceJALE(0:nalesteps, ie, 0:je, ke), &
672  flowdoms(nn, level, sps)%sFaceKALE(0:nalesteps, ie, je, 0:ke), stat=ierr)
673  if (ierr /= 0) &
674  call terminate("allocMemFlovarPart2", &
675  "Memory allocation failure for &
676  &sVeloIALE, sVeloJALE and sVeloKALE; &
677  &sFaceIALE, sFaceJALE and sFaceKALE.")
678  end if
679 
680  ! Test if we are on the finest mesh.
681 
682  fineleveltest: if (level == 1) then
683 
684  ! Allocate the memory that must always be allocated.
685 
686  allocate ( &
687  flowdoms(nn, level, sps)%dw(0:ib, 0:jb, 0:kb, 1:nw), &
688  flowdoms(nn, level, sps)%fw(0:ib, 0:jb, 0:kb, 1:nwf), &
689  flowdoms(nn, level, sps)%dtl(1:ie, 1:je, 1:ke), &
690  flowdoms(nn, level, sps)%radI(1:ie, 1:je, 1:ke), &
691  flowdoms(nn, level, sps)%radJ(1:ie, 1:je, 1:ke), &
692  flowdoms(nn, level, sps)%radK(1:ie, 1:je, 1:ke), &
693  flowdoms(nn, level, sps)%scratch(0:ib, 0:jb, 0:kb, 10), &
694  flowdoms(nn, level, sps)%shockSensor(0:ib, 0:jb, 0:kb), &
695  stat=ierr)
696  if (ierr /= 0) &
697  call terminate("allocMemFlovarPart2", &
698  "Memory allocation failure for dw, fw, dwOld, fwOld, &
699  &gamma, dtl and the spectral radii.")
700 
701  ! Initialize dw and fw to zero.
702 
703  flowdoms(nn, level, sps)%dw = zero
704  flowdoms(nn, level, sps)%fw = zero
705 
706  ! Extra variables for ALE
707  if (equationmode == unsteady .and. useale) then
708  allocate ( &
709  flowdoms(nn, level, sps)%dwALE(0:nalesteps, 0:ib, 0:jb, 0:kb, 1:nw), &
710  flowdoms(nn, level, sps)%fwALE(0:nalesteps, 0:ib, 0:jb, 0:kb, 1:nwf), &
711  stat=ierr)
712  if (ierr /= 0) &
713  call terminate("allocMemFlovarPart2", &
714  "Memory allocation failure for dwALE, fwALE.")
715 
716  flowdoms(nn, level, sps)%dwALE = zero
717  flowdoms(nn, level, sps)%fwALE = zero
718  end if
719 
720  ! Allocate the memory for the zeroth runge kutta stage
721  allocate (flowdoms(nn, level, sps)%wn(2:il, 2:jl, 2:kl, 1:nwf), &
722  flowdoms(nn, level, sps)%pn(2:il, 2:jl, 2:kl), stat=ierr)
723  if (ierr /= 0) &
724  call terminate("allocMemFlovarPart2", &
725  "Memory allocation failure for wn and pn")
726 
727  ! For unsteady mode using Runge-Kutta schemes allocate the
728  ! memory for dwOldRK.
729 
730  if (equationmode == unsteady .and. &
732 
733  mm = nrkstagesunsteady - 1
734  allocate (flowdoms(nn, level, sps)%dwOldRK(mm, il, jl, kl, nw), &
735  stat=ierr)
736  if (ierr /= 0) &
737  call terminate("allocMemFlovarPart2", &
738  "Memory allocation failure for dwOldRK.")
739  end if
740 
741  else fineleveltest
742 
743  ! Coarser level. Allocate the memory for the multigrid
744  ! forcing term and the state variables upon entry.
745 
746  allocate (flowdoms(nn, level, sps)%p1(1:ie, 1:je, 1:ke), &
747  flowdoms(nn, level, sps)%w1(1:ie, 1:je, 1:ke, 1:nwf), &
748  flowdoms(nn, level, sps)%wr(2:il, 2:jl, 2:kl, 1:nwf), &
749  stat=ierr)
750  if (ierr /= 0) &
751  call terminate("allocMemFlovarPart2", &
752  "Memory allocation failure for p1, w1 &
753  &and wr")
754 
755  ! Initialize w1 and p1 to zero, just that they
756  ! are initialized.
757 
758  flowdoms(nn, level, sps)%p1 = zero
759  flowdoms(nn, level, sps)%w1 = zero
760 
761  end if fineleveltest
762 
763  end do domains
764 
765  end subroutine allocmemflovarpart2
766 
767  subroutine allocrestartfiles(nFiles)
768  !
769  ! Allocate memory for the restartfles *
770  ! The array is populated from Python using setRestartFiles
771  ! If memeory has been allocated for the array there exist at
772  ! least one element in the array.
773  !
774  use constants
775  use inputio, only: restartfiles
776  use utils, only: terminate
777  implicit none
778  !
779  ! Subroutine argument.
780  !
781  integer(kind=intType) :: nFiles
782  !
783  ! Local variables.
784  !
785  integer :: ierr
786 
787  if (allocated(restartfiles)) then
788  deallocate (restartfiles)
789  end if
790 
791  allocate (restartfiles(nfiles), stat=ierr)
792  if (ierr /= 0) &
793  call terminate("allocRestartFiles", &
794  "Memory allocation failure for restartFiles")
795 
796  ! Zero Array with empty strings
797  restartfiles = ""
798 
799  end subroutine allocrestartfiles
800 
802  !
803  ! copySpectralSolution copies the solution of the 1st spectral
804  ! solution to all spectral solutions. This typically occurs when
805  ! a for the spectral mode a restart is made from a steady or an
806  ! unsteady solution. Possible rotation effects are taken into
807  ! account for the velocity components.
808  !
809  use constants
810  use block, only: flowdoms, ndom
811  use flowvarrefstate, only: nw
813  use iomodule, only: iovar
814  use monitor, only: timeunsteadyrestart
815  use section, only: sections, nsections
816  use utils, only: rotmatrixrigidbody
817  implicit none
818  !
819  ! Local variables.
820  !
821  integer(kind=intType) :: sps, spsm1, mm, nn, i, j, k, l
822 
823  real(kind=realtype) :: dt, tnew, told, tmp
824  real(kind=realtype) :: theta, costheta, sintheta
825 
826  real(kind=realtype), dimension(3) :: rotpoint, uu
827  real(kind=realtype), dimension(3) :: xt, yt, zt
828  real(kind=realtype), dimension(3, 3) :: rotmat
829 
830  real(kind=realtype), dimension(nSections, 3, 3) :: rotmatsec
831 
832  ! Determine the rotation matrix from one spectral solution to the
833  ! other for every section.
834 
835  sectionloop: do nn = 1, nsections
836 
837  ! Test if the section is rotating.
838 
839  testrotating: if (sections(nn)%rotating) then
840 
841  ! Section is rotating. Determine the angle between the
842  ! spectral solutions and its sine and cosine.
843 
844  i = sections(nn)%nSlices * ntimeintervalsspectral
845  theta = two * pi / real(i, realtype)
846  costheta = cos(theta)
847  sintheta = sin(theta)
848 
849  ! Transform to a frame where the xt-axis points in the
850  ! direction of the rotation vector.
851 
852  xt(1) = sections(nn)%rotAxis(1)
853  xt(2) = sections(nn)%rotAxis(2)
854  xt(3) = sections(nn)%rotAxis(3)
855 
856  ! Construct the yt axis. It does not matter exactly as long
857  ! as it is normal to xt.
858 
859  if (abs(xt(2)) < 0.707107_realtype) then
860  yt(1) = zero
861  yt(2) = one
862  yt(3) = zero
863  else
864  yt(1) = zero
865  yt(2) = zero
866  yt(3) = one
867  end if
868 
869  ! Make sure that yt is normal to xt.
870 
871  tmp = xt(1) * yt(1) + xt(2) * yt(2) + xt(3) * yt(3)
872  yt(1) = yt(1) - tmp * xt(1)
873  yt(2) = yt(2) - tmp * xt(2)
874  yt(3) = yt(3) - tmp * xt(3)
875 
876  ! And create a unit vector.
877 
878  tmp = one / sqrt(yt(1)**2 + yt(2)**2 + yt(3)**2)
879  yt(1) = tmp * yt(1)
880  yt(2) = tmp * yt(2)
881  yt(3) = tmp * yt(3)
882 
883  ! Create the vector zt by taking the cross product xt*yt.
884 
885  zt(1) = xt(2) * yt(3) - xt(3) * yt(2)
886  zt(2) = xt(3) * yt(1) - xt(1) * yt(3)
887  zt(3) = xt(1) * yt(2) - xt(2) * yt(1)
888 
889  ! The rotation matrix in the xt,yt,zt frame is given by
890  !
891  ! R = | 1 0 0 |
892  ! | 0 cos(theta) -sin(theta) |
893  ! | 0 sin(theta) cos(theta) |
894  !
895  ! The rotation matrix in the standard cartesian frame is then
896  ! given by t * r * t^t, where the colums of the transformation
897  ! matrix t are the unit vectors xt,yt,zt. One can easily check
898  ! this by checking rotation around the y- and z-axis. The
899  ! result of this is the expression below.
900 
901  rotmatsec(nn, 1, 1) = xt(1) * xt(1) &
902  + costheta * (yt(1) * yt(1) + zt(1) * zt(1))
903  rotmatsec(nn, 1, 2) = xt(1) * xt(2) &
904  + costheta * (yt(1) * yt(2) + zt(1) * zt(2)) &
905  - sintheta * (yt(1) * zt(2) - yt(2) * zt(1))
906  rotmatsec(nn, 1, 3) = xt(1) * xt(3) &
907  + costheta * (yt(1) * yt(3) + zt(1) * zt(3)) &
908  - sintheta * (yt(1) * zt(3) - yt(3) * zt(1))
909 
910  rotmatsec(nn, 2, 1) = xt(1) * xt(2) &
911  + costheta * (yt(1) * yt(2) + zt(1) * zt(2)) &
912  + sintheta * (yt(1) * zt(2) - yt(2) * zt(1))
913  rotmatsec(nn, 2, 2) = xt(2) * xt(2) &
914  + costheta * (yt(2) * yt(2) + zt(2) * zt(2))
915  rotmatsec(nn, 2, 3) = xt(2) * xt(3) &
916  + costheta * (yt(2) * yt(3) + zt(2) * zt(3)) &
917  - sintheta * (yt(2) * zt(3) - yt(3) * zt(2))
918 
919  rotmatsec(nn, 3, 1) = xt(1) * xt(3) &
920  + costheta * (yt(1) * yt(3) + zt(1) * zt(3)) &
921  + sintheta * (yt(1) * zt(3) - yt(3) * zt(1))
922  rotmatsec(nn, 3, 2) = xt(2) * xt(3) &
923  + costheta * (yt(2) * yt(3) + zt(2) * zt(3)) &
924  + sintheta * (yt(2) * zt(3) - yt(3) * zt(2))
925  rotmatsec(nn, 3, 3) = xt(3) * xt(3) &
926  + costheta * (yt(3) * yt(3) + zt(3) * zt(3))
927 
928  else testrotating
929 
930  ! Section is not rotating. Set the rotation matrix to the
931  ! identity matrix.
932 
933  rotmatsec(nn, 1, 1) = one
934  rotmatsec(nn, 1, 2) = zero
935  rotmatsec(nn, 1, 3) = zero
936 
937  rotmatsec(nn, 2, 1) = zero
938  rotmatsec(nn, 2, 2) = one
939  rotmatsec(nn, 2, 3) = zero
940 
941  rotmatsec(nn, 3, 1) = zero
942  rotmatsec(nn, 3, 2) = zero
943  rotmatsec(nn, 3, 3) = one
944 
945  end if testrotating
946 
947  end do sectionloop
948 
949  ! Initialize told to timeUnsteadyRestart. This takes the
950  ! possibility into account that the spectral mode is restarted
951  ! from an unsteady computation. Although not likely this
952  ! possibility is allowed and should therefore be taken into
953  ! account. Anyway told corresponds to the time of the 1st
954  ! spectral solution. Also determine the time step between the
955  ! spectral solutions. Both told and dt are only used to determine
956  ! the rigid body motion and if these are specified it is assumed
957  ! that there is only one section present in the grid.
958 
959  told = timeunsteadyrestart
960  dt = sections(1)%timePeriod &
961  / real(ntimeintervalsspectral, realtype)
962 
963  ! Loop over the number of spectral modes, starting at 2.
964 
965  spectralloop: do sps = 2, ntimeintervalsspectral
966 
967  ! Determine the corresponding time for this spectral solution
968  ! and store sps - 1 a bit easier.
969 
970  tnew = told + dt
971  spsm1 = sps - 1
972 
973  ! Determine the rotation matrix and rotation point between the
974  ! told and tnew for the rigid body rotation of the entire mesh.
975  ! The rotation point is not needed for the transformation of the
976  ! velocities, but rotMatrixRigidBody happens to compute it.
977 
978  call rotmatrixrigidbody(tnew, told, rotmat, rotpoint)
979 
980  ! Loop over the local number of blocks.
981 
982  domains: do nn = 1, ndom
983 
984  ! Store the section ID of this block a bit easier in mm.
985 
986  mm = flowdoms(nn, 1, 1)%sectionId
987 
988  ! Loop over the owned cells of this block. As the number of
989  ! cells is identical for all spectral solutions, it does not
990  ! matter which mode is taken for the upper dimensions.
991 
992  do k = 2, flowdoms(nn, 1, 1)%kl
993  do j = 2, flowdoms(nn, 1, 1)%jl
994  do i = 2, flowdoms(nn, 1, 1)%il
995 
996  ! Step 1. Copy the solution variables w from
997  ! the previous spectral solution.
998 
999  do l = 1, nw
1000  iovar(nn, sps)%w(i, j, k, l) = iovar(nn, spsm1)%w(i, j, k, l)
1001  end do
1002 
1003  ! Step 2. Apply the rigid body motion rotation matrix
1004  ! to the velocity. Use uu as a temporary storage.
1005 
1006  uu(1) = rotmat(1, 1) * iovar(nn, sps)%w(i, j, k, ivx) &
1007  + rotmat(1, 2) * iovar(nn, sps)%w(i, j, k, ivy) &
1008  + rotmat(1, 3) * iovar(nn, sps)%w(i, j, k, ivz)
1009 
1010  uu(2) = rotmat(2, 1) * iovar(nn, sps)%w(i, j, k, ivx) &
1011  + rotmat(2, 2) * iovar(nn, sps)%w(i, j, k, ivy) &
1012  + rotmat(2, 3) * iovar(nn, sps)%w(i, j, k, ivz)
1013 
1014  uu(3) = rotmat(3, 1) * iovar(nn, sps)%w(i, j, k, ivx) &
1015  + rotmat(3, 2) * iovar(nn, sps)%w(i, j, k, ivy) &
1016  + rotmat(3, 3) * iovar(nn, sps)%w(i, j, k, ivz)
1017 
1018  ! Step 3. Apply the rotation matrix of the section to
1019  ! the velocity.
1020 
1021  iovar(nn, sps)%w(i, j, k, ivx) = rotmatsec(mm, 1, 1) * uu(1) &
1022  + rotmatsec(mm, 1, 2) * uu(2) &
1023  + rotmatsec(mm, 1, 3) * uu(3)
1024 
1025  iovar(nn, sps)%w(i, j, k, ivy) = rotmatsec(mm, 2, 1) * uu(1) &
1026  + rotmatsec(mm, 2, 2) * uu(2) &
1027  + rotmatsec(mm, 2, 3) * uu(3)
1028 
1029  iovar(nn, sps)%w(i, j, k, ivz) = rotmatsec(mm, 3, 1) * uu(1) &
1030  + rotmatsec(mm, 3, 2) * uu(2) &
1031  + rotmatsec(mm, 3, 3) * uu(3)
1032  end do
1033  end do
1034  end do
1035 
1036  end do domains
1037 
1038  ! Set told to tnew for the next spectral solution.
1039 
1040  told = tnew
1041 
1042  end do spectralloop
1043 
1044  end subroutine copyspectralsolution
1045 
1047  !
1048  ! determineSolFileNames determines the number and names of the
1049  ! files that contain the solutions. For steady computations only
1050  ! one file must be present. For unsteady the situation is a
1051  ! little more complicated. It is attempted to read as many
1052  ! solutions as needed for a consistent restart. If not possible
1053  ! as many as possible solutions are read. For an unsteady
1054  ! computation the order will be reduced; for time spectral mode
1055  ! the solution will be interpolated.
1056  !
1057  use constants
1058  use communication, only: myid
1059  use inputio, only: restartfiles
1060  use inputphysics, only: equationmode
1064  use utils, only: terminate
1065  implicit none
1066  !
1067  ! Local variables
1068  !
1069  integer :: ierr
1070 
1071  integer(kind=intType) :: ii, nn
1072 
1073  character(len=7) :: integerString
1074  character(len=maxStringLen) :: tmpName
1075 
1076  ! Initialize copySpectral and interpolSpectral to .false.
1077 
1078  copyspectral = .false.
1079  interpolspectral = .false.
1080 
1081  ! Determine the desired number of files to be read. This depends
1082  ! on the equation mode we have to solve for. Also set the
1083  ! corresponding file names.
1084 
1085  select case (equationmode)
1086 
1087  case (steady)
1088 
1089  ! Steady computation. Only one solution needs to be read.
1090  ! In case a list of restart files were provided in the python
1091  ! script, we force a read from only the first solution file.
1092 
1093  nsolsread = 1
1094  allocate (solfiles(nsolsread), stat=ierr)
1095  if (ierr /= 0) &
1096  call terminate("determineSolFileNames", &
1097  "Memory allocation failure for solFiles")
1098 
1099  solfiles(1) = restartfiles(1)
1100 
1101  ! Check if the files can be opened, exit if that fails.
1102  call checksolfilenames()
1103 
1104  !===============================================================
1105 
1106  case (unsteady)
1107 
1108  ! Unsteady computation. For a consistent restart nOldLevels
1109  ! solutions must be read. All restart files are provided explicitly
1110  ! from python script.
1111  call setsolfilenames()
1112 
1113  ! Check if the files can be opened, exit if that fails.
1114  call checksolfilenames()
1115 
1116  ! Set nOldSolAvail to nSolsRead and check if a consistent
1117  ! restart can be made. If not, processor 0 prints a warning.
1118 
1120  if (myid == 0 .and. noldsolavail < noldlevels) then
1121 
1122  print "(a)", "#"
1123  print "(a)", "# Warning"
1124  print "(a)", "# Not enough data found for a consistent &
1125  &time accurate restart."
1126  print "(a)", "# Order is reduced in the first time steps &
1127  &until enough data is available again."
1128  print "(a)", "#"
1129 
1130  end if
1131 
1132  ! Set the logicals oldSolWritten, such that nothing is
1133  ! overwritten when solution files are dumped for this
1134  ! computation.
1135 
1136  ii = min(nsolsread, noldlevels - 1)
1137  do nn = 1, ii
1138  oldsolwritten(nn) = .true.
1139  end do
1140 
1141  !===============================================================
1142 
1143  case (timespectral)
1144 
1145  ! Time spectral computation. For a consistent restart
1146  ! nTimeIntervalsSpectral solutions must be read. First
1147  ! determine the the restart files.
1148  call setsolfilenames()
1149 
1150  ! Check if the files can be opened, exit if that fails.
1151  call checksolfilenames()
1152 
1153  ! Check whether or not the spectral solution must be copied
1154  ! or interpolated.
1155 
1156  if (nsolsread == 1) then
1157  copyspectral = .true.
1158  else if (nsolsread /= ntimeintervalsspectral) then
1159  interpolspectral = .true.
1160  end if
1161 
1162  end select
1163 
1164  end subroutine determinesolfilenames
1165 
1166  subroutine setsolfilenames
1167  !
1168  ! setSolFileNames allocates and set the solution files that
1169  ! will be read and loaded in the restart
1170  !
1171  use communication
1172  use inputio
1173  use utils, only: terminate
1174  use variablereading, only: solfiles, nsolsread
1175  implicit none
1176  !
1177  ! Local variables
1178  !
1179  integer :: ierr
1180  integer(kind=intType) :: nn
1181 
1182  ! The length of the array provided gives the number of nSolsRead.
1183  nsolsread = SIZE(restartfiles, 1)
1184 
1185  ! Allocate the memory for the file names and set them.
1186  allocate (solfiles(nsolsread), stat=ierr)
1187  if (ierr /= 0) &
1188  call terminate("determineSolFileNames", &
1189  "Memory allocation failure for solFiles")
1190 
1191  do nn = 1, nsolsread
1192  solfiles(nn) = restartfiles(nn)
1193  end do
1194 
1195  end subroutine setsolfilenames
1196 
1198  !
1199  ! checkSolFileNames will check if the provided restart files
1200  ! are readable on disk. If not readable return fail, if readable
1201  ! message will be printed to let the user know that restart file
1202  ! will be tried to read.
1203  !
1204  use communication
1205  use utils, only: terminate
1206  use variablereading, only: solfiles, nsolsread
1207  use commonformats, only: strings
1208  implicit none
1209  !
1210  ! Local variables
1211  !
1212  integer :: ierr
1213  character(len=maxStringLen) :: errorMessage
1214  integer(kind=intType) :: nn
1215 
1216  do nn = 1, nsolsread
1217  open (unit=21, file=solfiles(nn), status="old", iostat=ierr)
1218  if (ierr /= 0) then
1219  write (errormessage, *) "Restart file ", trim(solfiles(nn)), &
1220  " could not be opened for reading"
1221  call terminate("checkSolFileNames", errormessage)
1222  exit
1223  end if
1224  close (unit=21)
1225 
1226  if (myid == 0) then
1227  print "(a)", "#"
1228  write (*, strings) "# Found restart file: ", trim(solfiles(nn))
1229  print "(a)", "#"
1230  end if
1231  end do
1232 
1233  end subroutine checksolfilenames
1234 
1235  subroutine initdepvarandhalos(halosRead)
1236  !
1237  ! InitDepvarAndHalos computes the dependent flow variables,
1238  ! like viscosities, and initializes the halo cells by applying
1239  ! the boundary conditions and exchanging the internal halo's.
1240  ! This is all done on the start level grid.
1241  !
1242  use blockpointers
1243  use flowvarrefstate
1244  use inputio
1245  use inputiteration
1246  use inputphysics
1247  use inputtimespectral
1248  use iteration
1249  use monitor
1250  use section
1251  use utils, only: setpointers
1252  use haloexchange, only: whalo2
1253  use turbutils, only: computeeddyviscosity
1254  use turbbcroutines, only: applyallturbbc
1257  use flowutils, only: computelamviscosity
1258  use bcroutines, only: applyallbc
1259  use residuals, only: residual
1260  implicit none
1261  !
1262  ! Subroutine arguments.
1263  !
1264  logical, intent(in) :: halosRead
1265  !
1266  ! Local variables.
1267  !
1268  integer(kind=intType) :: nn, mm
1269  real(kind=realtype) :: relaxbleedsor
1270 
1271  real(kind=realtype), dimension(nSections) :: t
1272 
1273  logical :: initBleeds
1274 
1275  ! Set the logical whether or not to initialize the prescribed
1276  ! data for the bleed regions. If the halos were read the bleeds
1277  ! have been initialized already and nothing needs to be done.
1278 
1279  if (halosread) then
1280  initbleeds = .false.
1281  else
1282  initbleeds = .true.
1283  end if
1284 
1285  ! Compute the face velocities and for viscous walls the slip
1286  ! velocities. This is done for all the mesh levels.
1287 
1288  currentlevel = 1
1289  groundlevel = 1
1290 
1291  do mm = 1, ntimeintervalsspectral
1292 
1293  ! Compute the time, which corresponds to this spectral solution.
1294  ! For steady and unsteady mode this is simply the restart time;
1295  ! for the spectral mode the periodic time must be taken into
1296  ! account, which can be different for every section.
1297 
1299 
1300  if (equationmode == timespectral) then
1301  do nn = 1, nsections
1302  t(nn) = t(nn) + (mm - 1) * sections(nn)%timePeriod &
1303  / real(ntimeintervalsspectral, realtype)
1304  end do
1305  end if
1306  call gridvelocitiesfinelevel(.false., t, mm)
1308  call normalvelocitiesalllevels(mm)
1309 
1310  call slipvelocitiesfinelevel(.false., t, mm)
1312 
1313  end do
1314 
1315  ! Loop over the number of spectral solutions and blocks
1316  ! to compute the laminar viscosity.
1317 
1318  do mm = 1, ntimeintervalsspectral
1319  do nn = 1, ndom
1320  call setpointers(nn, mgstartlevel, mm)
1321  call computelamviscosity(.false.)
1322  end do
1323  end do
1324 
1325  ! Exchange the solution on the multigrid start level.
1326  ! It is possible that the halo values are needed for the boundary
1327  ! conditions. Viscosities are not exchanged.
1328 
1329  call whalo2(mgstartlevel, 1_inttype, nw, .true., .true., &
1330  .false.)
1331 
1332  ! Apply all flow boundary conditions to be sure that the halo's
1333  ! contain the correct values. These might be needed to compute
1334  ! the eddy-viscosity. Also the data for the outflow bleeds
1335  ! is determined.
1336 
1339 
1340  call applyallbc(.true.)
1341 
1342  ! Loop over the number of spectral solutions and blocks
1343  ! to compute the eddy viscosities.
1344 
1345  do mm = 1, ntimeintervalsspectral
1346  do nn = 1, ndom
1347 
1348  ! Set the pointers for this block.
1349 
1350  call setpointers(nn, mgstartlevel, mm)
1351 
1352  ! Compute the eddy viscosity for rans computations using
1353  ! an eddy viscosity model.
1354 
1355  call computeeddyviscosity(.false.)
1356 
1357  ! In case of a rans computation and no restart, initialize
1358  ! the turbulent variables a bit better for some turbulence
1359  ! models.
1360 
1361  ! if(equations == RANSEquations .and. (.not. restart)) then
1362  !
1363  ! select case (turbModel)
1364  !
1365  ! case (komegaWilcox, komegaModified, menterSST)
1366  ! call initKOmega(0_intType)
1367  ! call computeEddyViscosity
1368  !
1369  ! end select
1370  !
1371  ! endif
1372 
1373  end do
1374  end do
1375 
1376  ! Exchange the laminar and eddy viscosities.
1377 
1378  call whalo2(mgstartlevel, 1_inttype, 0_inttype, .false., &
1379  .false., .true.)
1380 
1381  if (equations == ransequations) then
1382  call applyallturbbc(.true.)
1383  end if
1384  call applyallbc(.true.)
1385 
1386  ! Exchange the solution for the second time to be sure that all
1387  ! halo's are initialized correctly. As this is the initialization
1388  ! phase, this is not critical.
1389 
1390  call whalo2(mgstartlevel, 1_inttype, nw, .true., .true., &
1391  .true.)
1392 
1393  end subroutine initdepvarandhalos
1394 
1395  subroutine initflowrestart
1396  !
1397  ! initFlowRestart loads restart information from the restart
1398  ! file into the state variables.
1399  !
1400  use constants
1401  use iomodule, only: iovar
1403  halosread
1404  use utils, only: terminate
1405  implicit none
1406  !
1407  ! Local variables.
1408  !
1409  integer(kind=intType) :: ierr
1410 
1411  ! Initialize halosRead to .false. This will be overwritten
1412  ! if halo values are read during the restart.
1413 
1414  halosread = .false.
1415 
1416  ! Determine the number and names of the solution files and
1417  ! determine the contents of IOVar, the data structure to
1418  ! generalize the IO.
1419 
1421  call setiovar
1422 
1423  ! Determine the format of the files and read them.
1424  ! Note that halosRead is possibly overwritten in the
1425  ! folloing select case statement below
1426 
1427  call readrestartfile()
1428 
1429  ! Copy or interpolate the spectral solution, if needed.
1430 
1433 
1434  ! Release the memory of the file names and IOVar.
1435 
1436  deallocate (solfiles, iovar, stat=ierr)
1437  if (ierr /= 0) &
1438  call terminate("initFlow", &
1439  "Deallocation failure for solFiles and IOVar")
1440 
1441  ! At the moment the pressure is stored at the location of the
1442  ! total energy. Copy the pressure to its own arrays and
1443  ! compute the total energy.
1444 
1446 
1447  ! Initialize the halo cells if a restart is performed and the
1448  ! entire flow field if this is not the case.
1449 
1451 
1452  ! Initialize the dependent flow variables and the halo values.
1454 
1455  end subroutine initflowrestart
1456 
1457  subroutine initflowfield
1458  !
1459  ! initFlowfield initializes the flow field to a uniform flow on
1460  ! the start level grid. Exception may be some turbulence
1461  ! variables, which are initialized a bit smarter.
1462  !
1463  use constants
1464  use communication, only: myid
1465  use inputiteration, only: ncycles, nsgstartup
1466  use inputphysics, only: equationmode
1467  use inputunsteady, only: ntimestepsfine
1468  use iteration, only: noldsolavail
1471  implicit none
1472  !
1473  ! Local variables.
1474  !
1475  integer(kind=intType) :: nn
1476 
1477  ! Initialize nTimeStepsRestart to 0 (no restart
1478  ! is performed) and allocate the memory for the arrays to store
1479  ! convergence history. This allocation is only to be done by
1480  ! processor 0. For an unsteady computation the entire convergence
1481  ! history is stored and the values after every time step is
1482  ! stored in a separate array.
1483 
1484  ntimestepsrestart = 0
1486 
1487  nn = nsgstartup + ncycles
1488  if (equationmode == unsteady) nn = ntimestepsfine * nn
1489  if (myid == 0) call allocconvarrays(nn)
1490  if (equationmode == unsteady .and. myid == 0) &
1492 
1493  ! Set nOldSolAvail to 1, to indicate that an unsteady
1494  ! computation should be started with a lower order
1495  ! time integration scheme, because not enough states
1496  ! in the past are available.
1497 
1498  noldsolavail = 1
1499 
1500  ! Initialize the flow field to uniform flow.
1501 
1502  call setuniformflow
1503 
1504  end subroutine initflowfield
1505 
1506  subroutine initializehalos(halosRead)
1507  !
1508  ! initializeHalos sets the flow variables in the halo cells
1509  ! using a constant extrapolation. If the halos are read only the
1510  ! second halos are initialized, otherwise both.
1511  !
1512  use constants
1513  use blockpointers, only: ndom, w, p, rlv, ib, jb, kb, &
1514  il, ie, jl, je, kl, ke, rev
1515  use flowvarrefstate, only: viscous, nw, eddymodel, muinf
1516  use inputiteration, only: mgstartlevel
1517  use inputphysics, only: eddyvisinfratio
1519  use utils, only: setpointers
1520  implicit none
1521  !
1522  ! Subroutine arguments.
1523  !
1524  logical, intent(in) :: halosRead
1525  !
1526  ! Local variables.
1527  !
1528  integer(kind=intType) :: nn, mm, i, j, k, l
1529  integer(kind=intType) :: jj, kk
1530 
1531  ! Loop over the number of spectral solutions and blocks.
1532 
1533  spectralloop: do mm = 1, ntimeintervalsspectral
1534  domains: do nn = 1, ndom
1535 
1536  ! Set the pointers for this block.
1537 
1538  call setpointers(nn, mgstartlevel, mm)
1539 
1540  ! Determine the situation we are dealing with.
1541 
1542  testhalosread: if (halosread) then
1543 
1544  ! The first layer of halo cells have been read. Initialize
1545  ! the second layer from the first layer.
1546 
1547  ! Halo cells in the i-direction.
1548 
1549  do k = 0, kb
1550  kk = max(1_inttype, min(k, ke))
1551  do j = 0, jb
1552  jj = max(1_inttype, min(j, je))
1553 
1554  do l = 1, nw
1555  w(0, j, k, l) = w(1, jj, kk, l)
1556  w(ib, j, k, l) = w(ie, jj, kk, l)
1557  end do
1558 
1559  p(0, j, k) = p(1, jj, kk)
1560  p(ib, j, k) = p(ie, jj, kk)
1561 
1562  end do
1563  end do
1564 
1565  ! Halo cells in j-direction. Note that the i-halo's have
1566  ! already been set.
1567 
1568  do k = 0, kb
1569  kk = max(1_inttype, min(k, ke))
1570  do i = 1, ie
1571 
1572  do l = 1, nw
1573  w(i, 0, k, l) = w(i, 1, kk, l)
1574  w(i, jb, k, l) = w(i, je, kk, l)
1575  end do
1576 
1577  p(i, 0, k) = p(i, 1, kk)
1578  p(i, jb, k) = p(i, je, kk)
1579 
1580  end do
1581  end do
1582 
1583  ! Halo cells in k-direction. Note that the halo's in both
1584  ! i and j direction have already been set.
1585 
1586  do j = 1, je
1587  do i = 1, ie
1588 
1589  do l = 1, nw
1590  w(i, j, 0, l) = w(i, j, 1, l)
1591  w(i, j, kb, l) = w(i, j, ke, l)
1592  end do
1593 
1594  p(i, j, 0) = p(i, j, 1)
1595  p(i, j, kb) = p(i, j, ke)
1596 
1597  end do
1598  end do
1599 
1600  else testhalosread
1601 
1602  ! No halo cells have been read. Initialize both layers
1603  ! using the internal value.
1604 
1605  ! Halo cells in the i-direction.
1606 
1607  do k = 0, kb
1608  kk = max(2_inttype, min(k, kl))
1609  do j = 0, jb
1610  jj = max(2_inttype, min(j, jl))
1611 
1612  do l = 1, nw
1613  w(0, j, k, l) = w(2, jj, kk, l)
1614  w(1, j, k, l) = w(2, jj, kk, l)
1615  w(ie, j, k, l) = w(il, jj, kk, l)
1616  w(ib, j, k, l) = w(il, jj, kk, l)
1617  end do
1618 
1619  p(0, j, k) = p(2, jj, kk)
1620  p(1, j, k) = p(2, jj, kk)
1621  p(ie, j, k) = p(il, jj, kk)
1622  p(ib, j, k) = p(il, jj, kk)
1623 
1624  end do
1625  end do
1626 
1627  ! Halo cells in j-direction. Note that the i-halo's have
1628  ! already been set.
1629 
1630  do k = 0, kb
1631  kk = max(2_inttype, min(k, kl))
1632  do i = 2, il
1633 
1634  do l = 1, nw
1635  w(i, 0, k, l) = w(i, 2, kk, l)
1636  w(i, 1, k, l) = w(i, 2, kk, l)
1637  w(i, je, k, l) = w(i, jl, kk, l)
1638  w(i, jb, k, l) = w(i, jl, kk, l)
1639  end do
1640 
1641  p(i, 0, k) = p(i, 2, kk)
1642  p(i, 1, k) = p(i, 2, kk)
1643  p(i, je, k) = p(i, jl, kk)
1644  p(i, jb, k) = p(i, jl, kk)
1645 
1646  end do
1647  end do
1648 
1649  ! Halo cells in k-direction. Note that the halo's in both
1650  ! i and j direction have already been set.
1651 
1652  do j = 2, jl
1653  do i = 2, il
1654 
1655  do l = 1, nw
1656  w(i, j, 0, l) = w(i, j, 2, l)
1657  w(i, j, 1, l) = w(i, j, 2, l)
1658  w(i, j, ke, l) = w(i, j, kl, l)
1659  w(i, j, kb, l) = w(i, j, kl, l)
1660  end do
1661 
1662  p(i, j, 0) = p(i, j, 2)
1663  p(i, j, 1) = p(i, j, 2)
1664  p(i, j, ke) = p(i, j, kl)
1665  p(i, j, kb) = p(i, j, kl)
1666 
1667  end do
1668  end do
1669 
1670  end if testhalosread
1671 
1672  ! Initialize the laminar and eddy viscosity, if appropriate,
1673  ! such that no uninitialized memory is present.
1674  ! As the viscosities are dependent variables their values
1675  ! are not read during a restart.
1676 
1677  if (viscous) rlv = muinf
1679 
1680  end do domains
1681  end do spectralloop
1682 
1683  end subroutine initializehalos
1685  !
1686  ! interpolateSpectralSolution uses a spectral interpolation to
1687  ! determine the initialization of the flow solution.
1688  ! The solution is interpolated from the solution read, which
1689  ! contains a different number of time instances and is stored in
1690  ! IOVar()%w. This variable can be found in IOModule.
1691  !
1692  use constants
1693  use blockpointers, only: w, il, jl, kl, ndom, sectionid
1694  use flowvarrefstate, only: nw
1696  use iomodule, only: iovar
1697  use section, only: sections, nsections
1699  use variablereading, only: nsolsread
1700  implicit none
1701  !
1702  ! Local variables.
1703  !
1704  integer :: ierr
1705 
1706  integer(kind=intType) :: jj, nn, ll, sps, i, j, k, l
1707 
1708  real(kind=realtype) :: t
1709 
1710  real(kind=realtype), dimension(nSolsRead) :: alpscal
1711  real(kind=realtype), dimension(nSections, nSolsRead, 3, 3) :: alpmat
1712 
1713  ! Loop over the number of spectral solutions to be interpolated.
1714 
1715  spectralloop: do sps = 1, ntimeintervalsspectral
1716 
1717  ! Determine the ratio of the time of this solution and the
1718  ! periodic time.
1719 
1720  nn = sps - 1
1721  t = real(nn, realtype) / real(ntimeintervalsspectral, realtype)
1722 
1723  ! Determine the interpolation coefficients for both the scalar
1724  ! and the vector quantities.
1725 
1726  call spectralinterpolcoef(nsolsread, t, alpscal, alpmat)
1727 
1728  ! Loop over the local number of blocks.
1729 
1730  domains: do nn = 1, ndom
1731 
1732  ! Set the pointers for this block to the finest grid level.
1733 
1734  call setpointers(nn, 1_inttype, sps)
1735 
1736  ! Loop over the number of variables to be interpolated.
1737 
1738  varloop: do l = 1, nw
1739 
1740  ! Check if this is a velocity variable.
1741 
1742  veltest: if (l == ivx .or. l == ivy .or. l == ivz) then
1743 
1744  ! Velocity variable. Set ll, which is the row in the
1745  ! matrix coefficients.
1746 
1747  select case (l)
1748  case (ivx)
1749  ll = 1
1750  case (ivy)
1751  ll = 2
1752  case (ivz)
1753  ll = 3
1754  end select
1755 
1756  ! Loop over the owned cells to interpolate the variable.
1757 
1758  do k = 2, kl
1759  do j = 2, jl
1760  do i = 2, il
1761 
1762  ! Initialization to zero and loop over the number of
1763  ! spectral solutions used in the interpolation and
1764  ! update the variable accordingly. Note that for the
1765  ! vector variables the matrix coefficients must be
1766  ! used; these matrices can be different for the
1767  ! sections present in the grid.
1768 
1769  w(i, j, k, l) = zero
1770  do jj = 1, nsolsread
1771  w(i, j, k, l) = w(i, j, k, l) &
1772  + alpmat(sectionid, jj, ll, 1) &
1773  * iovar(nn, jj)%w(i, j, k, ivx) &
1774  + alpmat(sectionid, jj, ll, 2) &
1775  * iovar(nn, jj)%w(i, j, k, ivy) &
1776  + alpmat(sectionid, jj, ll, 3) &
1777  * iovar(nn, jj)%w(i, j, k, ivz)
1778  end do
1779 
1780  end do
1781  end do
1782  end do
1783 
1784  else veltest
1785 
1786  ! Scalar variable.
1787  ! Loop over the owned cells to interpolate the variable.
1788 
1789  do k = 2, kl
1790  do j = 2, jl
1791  do i = 2, il
1792 
1793  ! Initialization to zero and loop over the number of
1794  ! spectral solutions used in the interpolation and
1795  ! update the variable accordingly.
1796 
1797  w(i, j, k, l) = zero
1798  do jj = 1, nsolsread
1799  w(i, j, k, l) = w(i, j, k, l) &
1800  + alpscal(jj) * iovar(nn, jj)%w(i, j, k, l)
1801  end do
1802 
1803  end do
1804  end do
1805  end do
1806 
1807  end if veltest
1808 
1809  end do varloop
1810 
1811  end do domains
1812 
1813  end do spectralloop
1814 
1815  ! Release the memory of w of IOVar.
1816 
1817  do sps = 1, nsolsread
1818  do nn = 1, ndom
1819  deallocate (iovar(nn, sps)%w, stat=ierr)
1820  if (ierr /= 0) &
1821  call terminate("interpolateSpectralSolution", &
1822  "Deallocation failure for w.")
1823  end do
1824  end do
1825 
1826  end subroutine interpolatespectralsolution
1827 
1829  !
1830  ! releaseExtraMemBCs releases the extra memory allocated in
1831  ! allocMemBcdata. This additional memory was allocated, such
1832  ! that alternative boundary condition treatments can be handled
1833  ! in setBCDataFineGrid.
1834  !
1835  use constants
1836  use blockpointers, only: flowdoms, ndom, bctype, bcdata, nbocos
1838  use utils, only: setpointers, terminate
1839  implicit none
1840  !
1841  ! Local variables.
1842  !
1843  integer :: ierr
1844 
1845  integer(kind=intType) :: mm, nn, sps, level, nLevels, ii
1846 
1847  ! Determine the number of multigrid levels.
1848 
1849  nlevels = ubound(flowdoms, 2)
1850 
1851  ! Loop over the number of multigrid level, spectral solutions
1852  ! and local blocks.
1853 
1854  levelloop: do level = 1, nlevels
1855  spectralloop: do sps = 1, ntimeintervalsspectral
1856  domainsloop: do nn = 1, ndom
1857 
1858  ! Have the pointers in blockPointers point to the
1859  ! current block to make everything more readable.
1860 
1861  call setpointers(nn, level, sps)
1862 
1863  ! Loop over the number of boundary subfaces for this block.
1864 
1865  bocoloop: do mm = 1, nbocos
1866 
1867  ! Determine the boundary condition we are having here.
1868 
1869  inflowtype:select case(bctype(mm))
1870 
1871  case (subsonicinflow)
1872 
1873  ! Subsonic inflow. Determine the boundary condition
1874  ! treatment and release the accordingly. Note that
1875  ! the boundary condition treatment of the finest mesh
1876  ! must be used, because this data is not available yet
1877  ! on the coarse grid.
1878 
1879  ii = &
1880  flowdoms(nn, 1, sps)%BCData(mm)%subsonicInletTreatment
1881 
1882  select case (ii)
1883 
1884  case (totalconditions)
1885 
1886  ! Total conditions used. Release the memory for
1887  ! the density and velocities and nullify their
1888  ! pointers.
1889 
1890  deallocate (bcdata(mm)%rho, bcdata(mm)%velx, &
1891  bcdata(mm)%vely, bcdata(mm)%velz, &
1892  stat=ierr)
1893  if (ierr /= 0) &
1894  call terminate("releaseExtraMemBCs", &
1895  "Deallocation failure for rho, &
1896  &velx, vely and velz")
1897 
1898  nullify (bcdata(mm)%rho)
1899  nullify (bcdata(mm)%velx)
1900  nullify (bcdata(mm)%vely)
1901  nullify (bcdata(mm)%velz)
1902  !===================================================
1903 
1904  case (massflow)
1905 
1906  ! Full velocity vector and density prescribed at
1907  ! inlet boundaries. Release the memory for the
1908  ! total conditions and flow directions and nullify
1909  ! their pointers.
1910 
1911  deallocate (bcdata(mm)%ptInlet, &
1912  bcdata(mm)%ttInlet, &
1913  bcdata(mm)%htInlet, &
1914  bcdata(mm)%flowXdirInlet, &
1915  bcdata(mm)%flowYdirInlet, &
1916  bcdata(mm)%flowZdirInlet, stat=ierr)
1917  if (ierr /= 0) &
1918  call terminate("releaseExtraMemBCs", &
1919  "Deallocation failure for the &
1920  &total conditions.")
1921 
1922  nullify (bcdata(mm)%ptInlet)
1923  nullify (bcdata(mm)%ttInlet)
1924  nullify (bcdata(mm)%htInlet)
1925  nullify (bcdata(mm)%flowXdirInlet)
1926  nullify (bcdata(mm)%flowYdirInlet)
1927  nullify (bcdata(mm)%flowZdirInlet)
1928 
1929  end select
1930 
1931  end select inflowtype
1932 
1933  end do bocoloop
1934  end do domainsloop
1935  end do spectralloop
1936  end do levelloop
1937 
1938  end subroutine releaseextramembcs
1939 
1940  subroutine setiovar
1941  !
1942  ! setIOVar allocates the memory for the derived data type IOVar,
1943  ! which is simplifies the reading. If an interpolation must be
1944  ! performed for the time spectral method also the solution of
1945  ! this IO type is allocated. For all other cases the pointers of
1946  ! IOVar are set the the appropriate entries of flowDoms, with
1947  ! possible offset due to the usage of pointers.
1948  !
1949  use constants
1950  use block, only: flowdoms, ndom
1951  use flowvarrefstate, only: nw
1952  use inputphysics, only: equationmode
1953  use iomodule, only: iovar
1954  use utils, only: terminate
1956  implicit none
1957  !
1958  ! Local variables.
1959  !
1960  integer :: ierr
1961 
1962  integer(kind=intType) :: nn, mm, il, jl, kl
1963 
1964  ! Allocate the memory for IOVar.
1965 
1966  allocate (iovar(ndom, nsolsread), stat=ierr)
1967  if (ierr /= 0) &
1968  call terminate("setIOVar", &
1969  "Memory allocation failure for solRead")
1970 
1971  ! Determine the equation mode we are solving and set the pointers
1972  ! of coorRead accordingly, or even allocate the memory, if needed.
1973 
1974  select case (equationmode)
1975 
1976  case (steady)
1977 
1978  ! Steady computation. Only one solution needs to be read.
1979  ! Loop over the number of blocks and set the pointers.
1980  ! No pointer offset is needed.
1981 
1982  do nn = 1, ndom
1983  iovar(nn, 1)%pointerOffset = 0
1984  iovar(nn, 1)%w => flowdoms(nn, 1, 1)%w(1:, 1:, 1:, :)
1985  end do
1986 
1987  !===============================================================
1988 
1989  case (unsteady)
1990 
1991  ! Unsteady computation. The first solution should be stored in
1992  ! w. For the others the pointers point to wOld. As the
1993  ! starting indices of wOld are 2, a pointer shift takes place
1994  ! here. I know this is a pain in the butt, but that's what
1995  ! we have to live with.
1996 
1997  do nn = 1, ndom
1998  iovar(nn, 1)%pointerOffset = 0
1999  iovar(nn, 1)%w => flowdoms(nn, 1, 1)%w(1:, 1:, 1:, :)
2000 
2001  do mm = 2, nsolsread
2002  iovar(nn, mm)%pointerOffset = -1
2003  iovar(nn, mm)%w => flowdoms(nn, 1, 1)%wOld(mm - 1, 2:, 2:, 2:, :)
2004  end do
2005  end do
2006 
2007  !===============================================================
2008 
2009  case (timespectral)
2010 
2011  ! Time spectral mode. A further check is required.
2012 
2013  testallocsolread: if (interpolspectral) then
2014 
2015  ! A restart is performed using a different number of time
2016  ! instances than the previous computation. Consequently the
2017  ! solutions will be interpolated later on. Hence some
2018  ! additional storage is required for the solutions and thus
2019  ! the w variables of IOVar are allocated. No halo data is
2020  ! needed here. No pointer offset either due to the explicit
2021  ! allocation.
2022 
2023  do nn = 1, ndom
2024  il = flowdoms(nn, 1, 1)%il
2025  jl = flowdoms(nn, 1, 1)%jl
2026  kl = flowdoms(nn, 1, 1)%kl
2027 
2028  do mm = 1, nsolsread
2029  iovar(nn, mm)%pointerOffset = 0
2030 
2031  allocate (iovar(nn, mm)%w(2:il, 2:jl, 2:kl, nw), stat=ierr)
2032  if (ierr /= 0) &
2033  call terminate("setIOVar", &
2034  "Memory allocation failure for w")
2035  end do
2036 
2037  end do
2038 
2039  else testallocsolread
2040 
2041  ! A restart is made using either 1 solution or the correct
2042  ! number of solution instances. In both cases simply set
2043  ! the pointers of IOVar. No pointer offset is needed.
2044 
2045  do nn = 1, ndom
2046  do mm = 1, nsolsread
2047  iovar(nn, mm)%pointerOffset = 0
2048  iovar(nn, mm)%w => flowdoms(nn, 1, mm)%w(1:, 1:, 1:, :)
2049  end do
2050  end do
2051 
2052  end if testallocsolread
2053 
2054  end select
2055 
2056  end subroutine setiovar
2057 
2058  subroutine setpressureandcomputeenergy(halosRead)
2059  !
2060  ! Due to the usage of the variable IOVar, which generalizes the
2061  ! IO and leads to reuse of code, currently the pressure is
2062  ! stored at the position of rhoE. In this routine that data is
2063  ! copied to the pressure array and the total energy is computed.
2064  ! Note that this routine is only called when a restart is done.
2065  !
2066  use constants
2067  use blockpointers, only: ndom, p, w, il, jl, kl
2068  use flowvarrefstate, only: kpresent
2070  use utils, only: setpointers
2071  use flowutils, only: computeetotblock
2072  implicit none
2073  !
2074  ! Subroutine arguments.
2075  !
2076  logical, intent(in) :: halosRead
2077  !
2078  ! Local variables.
2079  !
2080  integer(kind=intType) :: sps, nn, nHalo
2081  integer(kind=intType) :: i, j, k
2082  integer(kind=intType) :: iBeg, iEnd, jBeg, jEnd, kBeg, kEnd
2083 
2084  ! Set the value of nHalo, depending whether or not the halo cells
2085  ! have been read from the restart file.
2086 
2087  nhalo = 0
2088  if (halosread) nhalo = 1
2089 
2090  ! Loop over the number of time instances and the local blocks.
2091  ! As this routine is only called when a restart is performed,
2092  ! the MG start level is the finest level.
2093 
2094  do sps = 1, ntimeintervalsspectral
2095  do nn = 1, ndom
2096 
2097  ! Set the pointers to the correct block. As this routine is
2098  ! only called when a restart is performed, the MG start level
2099  ! is the finest level.
2100 
2101  call setpointers(nn, 1_inttype, sps)
2102 
2103  ! Determine the range for which the pressure must be computed.
2104 
2105  ibeg = 2 - nhalo; jbeg = 2 - nhalo; kbeg = 2 - nhalo
2106  iend = il + nhalo; jend = jl + nhalo; kend = kl + nhalo
2107 
2108  ! Copy the pressure for the required cells.
2109 
2110  do k = kbeg, kend
2111  do j = jbeg, jend
2112  do i = ibeg, iend
2113  p(i, j, k) = w(i, j, k, irhoe)
2114  end do
2115  end do
2116  end do
2117 
2118  ! Compute the total energy as well.
2119 
2120  call computeetotblock(ibeg, iend, jbeg, jend, kbeg, kend, kpresent)
2121  end do
2122  end do
2123 
2124  end subroutine setpressureandcomputeenergy
2125  subroutine setrestartfiles(fileName, i)
2126  !
2127  ! Populates the restartfiles
2128  ! The array is populated from Python using setRestartFiles
2129  !
2130  use constants
2131  use inputio, only: restartfiles
2132  implicit none
2133  !
2134  ! Subroutine argument.
2135  !
2136  character(len=*), intent(inout) :: fileName
2137  integer(kind=intType) :: i
2138 
2139  restartfiles(i) = filename
2140 
2141  end subroutine setrestartfiles
2142 
2143  subroutine setuniformflow
2144  !
2145  ! setUniformFlow set the flow variables of all local blocks on
2146  ! the start level to the uniform flow field.
2147  !
2148  use constants
2149  use blockpointers, only: w, dw, fw, flowdoms, ib, jb, kb, &
2150  rev, rlv, ndom, bcdata, p
2151  use communication
2152  use flowvarrefstate, only: eddymodel, viscous, muinf, nw, nwf, &
2153  pinfcorr, winf
2154  use inputiteration, only: mgstartlevel
2157  use utils, only: setpointers
2158  implicit none
2159  !
2160  ! Local variables.
2161  !
2162  integer :: ierr
2163 
2164  integer(kind=intType) :: nn, mm, i, j, k, l
2165 
2166  real(kind=realtype) :: tmp
2167 
2168  real(kind=realtype), dimension(3) :: dirloc, dirglob
2169 
2170  ! Loop over the number of spectral solutions and blocks.
2171  spectralloop: do mm = 1, ntimeintervalsspectral
2172  domains: do nn = 1, ndom
2173 
2174  ! Set the pointers for this block.
2175 
2176  call setpointers(nn, mgstartlevel, mm)
2177 
2178  ! Set the w-variables to the ones of the uniform flow field.
2179  do l = 1, nw
2180  do k = 0, kb
2181  do j = 0, jb
2182  do i = 0, ib
2183  w(i, j, k, l) = winf(l)
2184  dw(i, j, k, l) = zero
2185  end do
2186  end do
2187  end do
2188  end do
2189  !set this here for a reinitialize flow to eliminate possible NAN's
2190  do l = 1, nwf
2191  do k = 0, kb
2192  do j = 0, jb
2193  do i = 0, ib
2194  fw(i, j, k, l) = zero
2195  end do
2196  end do
2197  end do
2198  end do
2199 
2200  ! Set the pressure.
2201 
2202  p = pinfcorr
2203 
2204  ! Initialize the laminar and eddy viscosity, if appropriate,
2205  ! such that no uninitialized memory is present.
2206 
2207  if (viscous) rlv = muinf
2209 
2210  end do domains
2211  end do spectralloop
2212 
2213  ! Correct for the time spectral method in combination with an
2214  ! internal flow computation the velocity direction.
2215  ! It is possible that the prescribed direction is different
2216  ! for every time instance.
2217 
2218  testcorrection: if (equationmode == timespectral .and. &
2219  flowtype == internalflow) then
2220 
2221  ! Loop over the number of spectral solutions.
2222 
2223  spectralloopcorr: do mm = 1, ntimeintervalsspectral
2224 
2225  ! Initialize the local direction to zero. In the loop
2226  ! below this direction will be accumulated.
2227 
2228  dirloc = zero
2229 
2230  ! Loop over the number of blocks and determine the average
2231  ! velocity direction prescribed at the inlets.
2232 
2233  domainloop1: do nn = 1, ndom
2234 
2235  ! Set the pointer for the BCData to make the code more
2236  ! readable.
2237 
2238  bcdata => flowdoms(nn, mgstartlevel, mm)%BCData
2239 
2240  ! Loop over the number of boundary subfaces and update
2241  ! the local prescribed velocity direction.
2242 
2243  do l = 1, flowdoms(nn, mgstartlevel, mm)%nBocos
2244  call velmagnanddirectionsubface(tmp, dirloc, bcdata, l)
2245  end do
2246  end do domainloop1
2247 
2248  ! Determine the sum of dirLoc and create a unit vector
2249  ! for the global direction.
2250 
2251  call mpi_allreduce(dirloc, dirglob, 3, adflow_real, mpi_sum, &
2252  adflow_comm_world, ierr)
2253 
2254  tmp = one / max(eps, sqrt(dirglob(1)**2 &
2255  + dirglob(2)**2 &
2256  + dirglob(3)**2))
2257  dirglob(1) = tmp * dirglob(1)
2258  dirglob(2) = tmp * dirglob(2)
2259  dirglob(3) = tmp * dirglob(3)
2260 
2261  ! Loop again over the local domains and correct the
2262  ! velocity direction.
2263 
2264  domainsloop2: do nn = 1, ndom
2265 
2266  ! Set the pointers for this block.
2267 
2268  call setpointers(nn, mgstartlevel, mm)
2269 
2270  do k = 0, kb
2271  do j = 0, jb
2272  do i = 0, ib
2273  tmp = sqrt(w(i, j, k, ivx)**2 &
2274  + w(i, j, k, ivy)**2 &
2275  + w(i, j, k, ivz)**2)
2276  w(i, j, k, ivx) = tmp * dirglob(1)
2277  w(i, j, k, ivy) = tmp * dirglob(2)
2278  w(i, j, k, ivz) = tmp * dirglob(3)
2279  end do
2280  end do
2281  end do
2282 
2283  end do domainsloop2
2284  end do spectralloopcorr
2285 
2286  end if testcorrection
2287 
2288  end subroutine setuniformflow
2289 
2290  !=================================================================
2291 
2292  subroutine velmagnanddirectionsubface(vmag, dir, BCData, mm)
2293  !
2294  ! VelMagnAndDirectionSubface determines the maximum value
2295  ! of the magnitude of the velocity as well as the sum of the
2296  ! flow directions for the currently active subface.
2297  !
2298  use constants
2299  use block
2300  implicit none
2301  !
2302  ! Subroutine arguments.
2303  !
2304  integer(kind=intType), intent(in) :: mm
2305 
2306  real(kind=realtype), intent(out) :: vmag
2307  real(kind=realtype), dimension(3), intent(inout) :: dir
2308 
2309  type(bcdatatype), dimension(:), pointer :: bcdata
2310  !
2311  ! Local variables.
2312  !
2313  integer(kind=intType) :: i, j
2314  real(kind=realtype) :: vel
2315 
2316  ! Initialize vmag to -1.0.
2317 
2318  vmag = -one
2319 
2320  ! Check if the velocity is prescribed.
2321 
2322  if (associated(bcdata(mm)%velx) .and. &
2323  associated(bcdata(mm)%vely) .and. &
2324  associated(bcdata(mm)%velz)) then
2325 
2326  ! Loop over the owned faces of the subface. As the cell range
2327  ! may contain halo values, the nodal range is used.
2328 
2329  do j = (bcdata(mm)%jnBeg + 1), bcdata(mm)%jnEnd
2330  do i = (bcdata(mm)%inBeg + 1), bcdata(mm)%inEnd
2331 
2332  ! Compute the magnitude of the velocity and compare it
2333  ! with the current maximum. Store the maximum of the two.
2334 
2335  vel = sqrt(bcdata(mm)%velx(i, j)**2 &
2336  + bcdata(mm)%vely(i, j)**2 &
2337  + bcdata(mm)%velz(i, j)**2)
2338  vmag = max(vmag, vel)
2339 
2340  ! Compute the unit vector of the velocity and add it to dir.
2341 
2342  vel = one / max(eps, vel)
2343  dir(1) = dir(1) + vel * bcdata(mm)%velx(i, j)
2344  dir(2) = dir(2) + vel * bcdata(mm)%vely(i, j)
2345  dir(3) = dir(3) + vel * bcdata(mm)%velz(i, j)
2346 
2347  end do
2348  end do
2349  end if
2350 
2351  ! Check if the velocity direction is prescribed.
2352 
2353  if (associated(bcdata(mm)%flowXdirInlet) .and. &
2354  associated(bcdata(mm)%flowYdirInlet) .and. &
2355  associated(bcdata(mm)%flowZdirInlet)) then
2356 
2357  ! Add the unit vectors to dir by looping over the owned
2358  ! faces of the subfaces. Again the nodal range must be
2359  ! used for this.
2360 
2361  do j = (bcdata(mm)%jnBeg + 1), bcdata(mm)%jnEnd
2362  do i = (bcdata(mm)%inBeg + 1), bcdata(mm)%inEnd
2363 
2364  dir(1) = dir(1) + bcdata(mm)%flowXdirInlet(i, j)
2365  dir(2) = dir(2) + bcdata(mm)%flowYdirInlet(i, j)
2366  dir(3) = dir(3) + bcdata(mm)%flowZdirInlet(i, j)
2367 
2368  end do
2369  end do
2370 
2371  end if
2372 
2373  end subroutine velmagnanddirectionsubface
2374  subroutine timespectralcoef(coefSpectral, matrixCoefSpectral, &
2375  diagMatCoefSpectral)
2376  !
2377  ! timeSpectralCoef computes the time integration coefficients
2378  ! for the time spectral method. As it is possible that sections
2379  ! have different periodic times these coefficients are
2380  ! determined for all the sections. For vector quantities, such
2381  ! as momentum, these coefficients can also be different due to
2382  ! rotation and the fact that only a part of the wheel is
2383  ! simulated.
2384  !
2385  use constants
2386  use flowvarrefstate, only: timeref
2388  use section, only: nsections, sections
2389  implicit none
2390  !
2391  ! Subroutine arguments.
2392  !
2393  real(kind=realtype), &
2394  dimension(nSections, nTimeIntervalsSpectral - 1), &
2395  intent(out) :: coefspectral
2396  real(kind=realtype), &
2397  dimension(nSections, nTimeIntervalsSpectral - 1, 3, 3), &
2398  intent(out) :: matrixcoefspectral
2399  real(kind=realtype), dimension(nSections, 3, 3), &
2400  intent(out) :: diagmatcoefspectral
2401  !
2402  ! Local variables.
2403  !
2404  integer(kind=intType) :: pp, nn, mm, ii, i, j, ntot
2405  real(kind=realtype) :: coef, dangle, angle, fact, slicesfact
2406 
2407  real(kind=realtype), dimension(3, 3) :: rotmat, tmp
2408 
2409  ! Loop over the number of sections.
2410 
2411  sectionloop: do mm = 1, nsections
2412 
2413  ! Initialize dAngle (smallest angle in the cotangent function)
2414  ! and coef, which is the multiplication factor in front of the
2415  ! cotangent/cosecant function. Coef is a combination of the 1/2
2416  ! and the 2*pi/timePeriod
2417 
2418  dangle = pi / real(ntimeintervalsspectral, realtype)
2419  coef = pi * timeref / sections(mm)%timePeriod
2420 
2421  ! Computation of the scalar coefficients.
2422 
2423  scalarloop: do nn = 1, (ntimeintervalsspectral - 1)
2424 
2425  angle = nn * dangle
2426 
2427  ! The coefficient for an odd and even number of time
2428  ! instances are different; the former is 1/sin, the
2429  ! latter cos/sin or 1/tan.
2430 
2431  coefspectral(mm, nn) = coef / sin(angle)
2432 
2433  if (mod(ntimeintervalsspectral, 2_inttype) == 0) &
2434  coefspectral(mm, nn) = coefspectral(mm, nn) * cos(angle)
2435 
2436  ! Negate coef for the next spectral coefficient.
2437 
2438  coef = -coef
2439 
2440  end do scalarloop
2441 
2442  ! Initialize dAngle to the smallest angle in the cotangent
2443  ! or cosecant function. Now this angle is for the entire wheel,
2444  ! i.e. the number of slices must be taken into account.
2445 
2446  ntot = ntimeintervalsspectral * sections(mm)%nSlices
2447  dangle = pi / real(ntot, realtype)
2448 
2449  ! Initialize the rotation matrix to the unity matrix.
2450 
2451  rotmat(1, 1) = one; rotmat(1, 2) = zero; rotmat(1, 3) = zero
2452  rotmat(2, 1) = zero; rotmat(2, 2) = one; rotmat(2, 3) = zero
2453  rotmat(3, 1) = zero; rotmat(3, 2) = zero; rotmat(3, 3) = one
2454 
2455  ! Loop over the number of spectral coefficient to initialize the
2456  ! matrix coefficients; this is basically pp == 0 in the loop
2457  ! over the number slices. Use is made of the fact that the
2458  ! rotation matrix is the identity for pp == 0.
2459  ! coef changes sign at every time instance
2460 
2461  slicesfact = one / real(sections(mm)%nSlices, realtype)
2462  fact = one
2463 
2464  do nn = 1, (ntimeintervalsspectral - 1)
2465 
2466  ! Determine the scalar coefficient. This value depends now
2467  ! whether the total number of time instances in the wheel is
2468  ! odd or even.
2469 
2470  angle = nn * dangle
2471  coef = one / sin(angle)
2472 
2473  if (mod(ntot, 2_inttype) == 0) &
2474  coef = coef * cos(angle)
2475 
2476  coef = coef * fact * slicesfact
2477 
2478  ! The first part of matrixCoefSpectral is a diagonal matrix,
2479  ! because this indicates the contribution of the current
2480  ! slice to the time derivative.
2481 
2482  matrixcoefspectral(mm, nn, 1, 1) = coef
2483  matrixcoefspectral(mm, nn, 1, 2) = zero
2484  matrixcoefspectral(mm, nn, 1, 3) = zero
2485 
2486  matrixcoefspectral(mm, nn, 2, 1) = zero
2487  matrixcoefspectral(mm, nn, 2, 2) = coef
2488  matrixcoefspectral(mm, nn, 2, 3) = zero
2489 
2490  matrixcoefspectral(mm, nn, 3, 1) = zero
2491  matrixcoefspectral(mm, nn, 3, 2) = zero
2492  matrixcoefspectral(mm, nn, 3, 3) = coef
2493 
2494  fact = -fact
2495 
2496  end do
2497 
2498  ! Initialize diagMatCoefSpectral to zero, because the
2499  ! starting index in the loop over the number of slices -1 is
2500  ! 1, i.e. the slice where the actual computation takes places
2501  ! does not contribute to diagMatCoefSpectral.
2502 
2503  do j = 1, 3
2504  do i = 1, 3
2505  diagmatcoefspectral(mm, i, j) = zero
2506  end do
2507  end do
2508 
2509  ! Loop over the additional slices which complete an entire
2510  ! revolution. To be able to compute the coefficients a bit
2511  ! easier the loop runs from 1 to nSlices-1 and not from
2512  ! 2 to nSlices.
2513 
2514  slicesloop: do pp = 1, (sections(mm)%nSlices - 1)
2515 
2516  ! Compute the rotation matrix for this slice. This is the
2517  ! old one multiplied by the transformation matrix going from
2518  ! one slices to the next. Use tmp as temporary storage.
2519 
2520  do j = 1, 3
2521  do i = 1, 3
2522  tmp(i, j) = rotmatrixspectral(mm, i, 1) * rotmat(1, j) &
2523  + rotmatrixspectral(mm, i, 2) * rotmat(2, j) &
2524  + rotmatrixspectral(mm, i, 3) * rotmat(3, j)
2525  end do
2526  end do
2527 
2528  rotmat = tmp
2529 
2530  slicesfact = one / real(sections(mm)%nSlices, realtype)
2531 
2532  ! Loop over the number of spectral coefficients and update
2533  ! matrixCoefSpectral. The multiplication with (-1)**nn
2534  ! takes place here too.
2535 
2536  ! Multiply also by the term (-1)**(pN+1)
2537 
2538  fact = one
2539  if (mod(pp * ntimeintervalsspectral, 2_inttype) /= 0) &
2540  fact = -one
2541  slicesfact = fact * slicesfact
2542 
2543  fact = one
2544  ii = pp * ntimeintervalsspectral
2545  do nn = 1, (ntimeintervalsspectral - 1)
2546 
2547  ! Compute the coefficient multiplying the rotation matrix.
2548  ! Again make a distinction between an odd and an even
2549  ! number of time instances for the entire wheel.
2550 
2551  angle = (nn + ii) * dangle
2552  coef = one / sin(angle)
2553 
2554  if (mod(ntot, 2_inttype) == 0) &
2555  coef = coef * cos(angle)
2556 
2557  coef = coef * fact * slicesfact
2558 
2559  ! Update matrixCoefSpectral.
2560 
2561  do j = 1, 3
2562  do i = 1, 3
2563  matrixcoefspectral(mm, nn, i, j) = &
2564  matrixcoefspectral(mm, nn, i, j) + coef * rotmat(i, j)
2565  end do
2566  end do
2567 
2568  fact = -fact
2569 
2570  end do
2571 
2572  ! Update diagMatCoefSpectral. Also here the distinction
2573  ! between odd and even number of time instances.
2574 
2575  angle = ii * dangle
2576  coef = one / sin(angle)
2577 
2578  if (mod(ntot, 2_inttype) == 0) &
2579  coef = coef * cos(angle)
2580 
2581  coef = coef * slicesfact
2582 
2583  do j = 1, 3
2584  do i = 1, 3
2585  diagmatcoefspectral(mm, i, j) = &
2586  diagmatcoefspectral(mm, i, j) - coef * rotmat(i, j)
2587  end do
2588  end do
2589 
2590  end do slicesloop
2591 
2592  ! The matrix coefficients must be multiplied by the leading
2593  ! coefficient, which depends on the actual periodic time.
2594 
2595  coef = pi * timeref / sections(mm)%timePeriod
2596 
2597  do j = 1, 3
2598  do i = 1, 3
2599  diagmatcoefspectral(mm, i, j) = &
2600  coef * diagmatcoefspectral(mm, i, j)
2601  end do
2602  end do
2603 
2604  do nn = 1, (ntimeintervalsspectral - 1)
2605 
2606  do j = 1, 3
2607  do i = 1, 3
2608  matrixcoefspectral(mm, nn, i, j) = &
2609  coef * matrixcoefspectral(mm, nn, i, j)
2610  end do
2611  end do
2612 
2613  end do
2614 
2615  end do sectionloop
2616 
2617  end subroutine timespectralcoef
2618 
2620  !
2621  ! timeSpectralMatrices computes the matrices for the time
2622  ! derivative of the time spectral method for all sections. For
2623  ! scalar quantities these matrices only differ if sections have
2624  ! different periodic times. For vector quantities, such as
2625  ! momentum, these matrices can be different depending on whether
2626  ! the section is rotating or not and the number of slices
2627  ! present.
2628  !
2629  use constants
2630  use inputphysics, only: equationmode
2633  use section, only: nsections
2634  use utils, only: terminate
2635  implicit none
2636  !
2637  ! Local variables.
2638  !
2639  integer :: ierr
2640 
2641  integer(kind=intType) :: nn, mm, ll, kk, ii
2642  integer(kind=intType) :: i, j
2643 
2644  real(kind=realtype), dimension(3, 3) :: tmpmat
2645 
2646  real(kind=realtype), dimension(:, :), allocatable :: coefspectral
2647  real(kind=realtype), dimension(:, :, :, :), allocatable :: &
2648  matrixcoefspectral
2649  real(kind=realtype), dimension(:, :, :), allocatable :: &
2650  diagmatcoefspectral
2651  !
2652  ! This routine is only used for the spectral solutions. Return
2653  ! immediately if a different mode is solved.
2654 
2655  if (equationmode /= timespectral) return
2656 
2657  ! Allocate the memory for the matrices as well as the help
2658  ! variables needed to construct these matrices.
2659 
2660  !added to allow second call in mdUpdateRoutines
2661  if (allocated(dscalar)) deallocate (dscalar)
2662  if (allocated(dvector)) deallocate (dvector)
2663 
2665  mm = 3 * nn
2666  kk = nn - 1
2667 
2668  allocate (dscalar(nsections, nn, nn), &
2669  dvector(nsections, mm, mm), &
2670  coefspectral(nsections, kk), &
2671  matrixcoefspectral(nsections, kk, 3, 3), &
2672  diagmatcoefspectral(nsections, 3, 3), stat=ierr)
2673  if (ierr /= 0) &
2674  call terminate("timeSpectralMatrices", &
2675  "Memory allocation failure for the matrices of &
2676  &the spectral time derivatives.")
2677 
2678  ! Determine the help variables needed to construct the
2679  ! actual matrices.
2680 
2681  call timespectralcoef(coefspectral, matrixcoefspectral, &
2682  diagmatcoefspectral)
2683  !
2684  ! Determine the time derivative matrices for the sections.
2685  !
2686  ! Loop over the number of sections.
2687 
2688  sectionloop: do ii = 1, nsections
2689  !
2690  ! Matrix for scalar quantities.
2691  !
2692  ! Loop over the number of rows.
2693 
2694  do nn = 1, ntimeintervalsspectral
2695 
2696  ! Set the diagonal element to zero, i.e. there is no
2697  ! contribution to the own time derivative.
2698 
2699  dscalar(ii, nn, nn) = zero
2700 
2701  ! Loop over the rest of the columns.
2702 
2703  do mm = 1, (ntimeintervalsspectral - 1)
2704 
2705  ! Determine the corresponding column index.
2706 
2707  ll = nn + mm
2708  if (ll > ntimeintervalsspectral) &
2709  ll = ll - ntimeintervalsspectral
2710 
2711  ! Store the corresponding coefficient in dscalar.
2712 
2713  dscalar(ii, nn, ll) = coefspectral(ii, mm)
2714 
2715  end do
2716  end do
2717  !
2718  ! Matrices for vector quantities.
2719  !
2720  ! Loop over the number of time intervals; the number of rows
2721  ! is 3 times this number.
2722 
2723  rowloop: do nn = 1, ntimeintervalsspectral
2724 
2725  ! Initialize the diagonal block to diagMatCoefSpectral,
2726  ! the additional diagonal entry needed for the rotational
2727  ! periodicity.
2728 
2729  kk = 3 * (nn - 1)
2730  do j = 1, 3
2731  do i = 1, 3
2732  dvector(ii, kk + i, kk + j) = diagmatcoefspectral(ii, i, j)
2733  end do
2734  end do
2735 
2736  ! Loop over the other time intervals, which contribute to
2737  ! the time derivative.
2738 
2739  columnloop: do mm = 1, (ntimeintervalsspectral - 1)
2740 
2741  ! Determine the corresponding column index and check the
2742  ! situation we are having here.
2743 
2744  ll = nn + mm
2745  if (ll > ntimeintervalsspectral) then
2746 
2747  ! Index is outside the range and a shift must be applied.
2748 
2749  ll = ll - ntimeintervalsspectral
2750 
2751  ! The vector must be rotated. This effect is incorporated
2752  ! directly in the matrix of time derivatives.
2753 
2754  do j = 1, 3
2755  do i = 1, 3
2756  tmpmat(i, j) = matrixcoefspectral(ii, mm, i, 1) &
2757  * rotmatrixspectral(ii, 1, j) &
2758  + matrixcoefspectral(ii, mm, i, 2) &
2759  * rotmatrixspectral(ii, 2, j) &
2760  + matrixcoefspectral(ii, mm, i, 3) &
2761  * rotmatrixspectral(ii, 3, j)
2762  end do
2763  end do
2764 
2765  else
2766 
2767  ! Index is in the range. Copy the matrix coefficient
2768  ! into tmpMat.
2769 
2770  do j = 1, 3
2771  do i = 1, 3
2772  tmpmat(i, j) = matrixcoefspectral(ii, mm, i, j)
2773  end do
2774  end do
2775 
2776  end if
2777 
2778  ! Determine the offset for the column index and store
2779  ! this submatrix in the correct place of dvector.
2780 
2781  ll = 3 * (ll - 1)
2782  do j = 1, 3
2783  do i = 1, 3
2784  dvector(ii, kk + i, ll + j) = tmpmat(i, j)
2785  end do
2786  end do
2787 
2788  end do columnloop
2789  end do rowloop
2790  end do sectionloop
2791 
2792  ! Release the memory of the help variables needed to construct
2793  ! the matrices of the time derivatives.
2794 
2795  deallocate (coefspectral, matrixcoefspectral, &
2796  diagmatcoefspectral, stat=ierr)
2797  if (ierr /= 0) &
2798  call terminate("timeSpectralMatrices", &
2799  "Deallocation failure for the help variables.")
2800 
2801  end subroutine timespectralmatrices
2802 
2803  subroutine readrestartfile()
2804  !
2805  ! readRestartFile reads the fine grid solution(s) from the
2806  ! restart file(s). If the restart file(s) do not correspond to
2807  ! the current mesh, the solution(s) are interpolated onto this
2808  ! mesh. It is also allowed to change boundary conditions, e.g.
2809  ! an alpha and/or Mach sweep is possible. Furthermore there is
2810  ! some support when starting from a different turbulence model,
2811  ! although this should be used with care.
2812  !
2813  use constants
2814  use cgnsgrid
2815  use su_cgns
2816  use variablereading ! Full import since we need basically everything
2817  use blockpointers, only: ibegor, jbegor, kbegor, il, jl, kl, ndom, &
2818  nbkglobal, nx, ny, nz
2820  use inputphysics, only: equationmode
2823  use utils, only: terminate, setpointers
2824  use sorting, only: bsearchstrings
2825  use commonformats, only: strings, stringint1
2826  implicit none
2827  !
2828  ! Local variables.
2829  !
2830  integer :: nzones, celldim, physdim, ierr, nsols
2831 
2832  integer(cgsize_t), dimension(9) :: sizes
2833  integer, dimension(9) :: rindsizes
2834  integer, dimension(nSolsRead) :: fileids
2835 
2836  integer(kind=intType) :: ii, jj, nn
2837  integer(kind=intType) :: ntypemismatch
2838  integer(kind=intType) :: nhimin, nhjmin, nhkmin
2839  integer(kind=intType) :: nhimax, nhjmax, nhkmax
2840 
2841  character(len=7) :: integerstring
2842  character(len=maxCGNSNameLen) :: cgnsname
2843  character(len=2*maxStringLen) :: errormessage
2844 
2845  ! Initialize halosRead to .true. This will be overwritten if
2846  ! there is at least one block present for which the halo data
2847  ! cannot be read.
2848 
2849  halosread = .true.
2850 
2851  ! Initialize nTypeMismatch to 0.
2852 
2853  ntypemismatch = 0
2854 
2855  ! Loop over the number of files to be read and open them.
2856 
2857  fileopenloop: do solid = 1, nsolsread
2858 
2859  ! Open the restart file for reading.
2860 
2861  call cg_open_f(solfiles(solid), mode_read, cgnsind, ierr)
2862  if (ierr /= all_ok) then
2863  write (errormessage, *) "File ", trim(solfiles(solid)), &
2864  " could not be opened for reading"
2865  call terminate("readRestartFile", errormessage)
2866  end if
2867 
2868  fileids(solid) = cgnsind
2869 
2870  ! Determine the number of bases in the cgns file.
2871  ! This must be at least 1.
2872 
2873  call cg_nbases_f(cgnsind, cgnsbase, ierr)
2874  if (ierr /= all_ok) &
2875  call terminate("readRestartFile", &
2876  "Something wrong when calling cg_nbases_f")
2877 
2878  if (cgnsbase < 1) then
2879  write (errormessage, *) "CGNS file ", trim(solfiles(solid)), &
2880  " does not contain a base"
2881  call terminate("readRestartFile", errormessage)
2882  end if
2883 
2884  ! Only data from the first base is read. Information from
2885  ! higher bases is ignored.
2886 
2887  cgnsbase = 1
2888 
2889  ! Read the cell and physical dimensions as well as the name for
2890  ! this base.
2891 
2892  call cg_base_read_f(cgnsind, cgnsbase, cgnsname, celldim, &
2893  physdim, ierr)
2894  if (ierr /= all_ok) &
2895  call terminate("readRestartFile", &
2896  "Something wrong when calling cg_base_read_f")
2897 
2898  ! Check the cell and physical dimensions. Both must be 3 for
2899  ! this code to work.
2900 
2901  if (celldim /= 3 .or. physdim /= 3) then
2902  write (errormessage, stringint1) "Both the number of cell and physical dimensions should be 3, not ", &
2903  celldim, " and ", physdim
2904  call terminate("readRestartFile", errormessage)
2905  end if
2906 
2907  end do fileopenloop
2908 
2909  ! Broadcast nTimeStepsRestart and timeUnsteadyRestart to all
2910  ! processors. These values are needed to perform a consistent
2911  ! unsteady restart.
2912 
2913  call mpi_bcast(ntimestepsrestart, 1, adflow_integer, 0, &
2914  adflow_comm_world, ierr)
2915  call mpi_bcast(timeunsteadyrestart, 1, adflow_real, 0, &
2916  adflow_comm_world, ierr)
2917 
2918  ! Get the scaling factors for density, pressure and velocity
2919  ! by reading the reference state.
2920 
2921  call scalefactors(fileids)
2922 
2923  ! Loop over the number of files to be read and read the solution.
2924 
2925  solloop: do solid = 1, nsolsread
2926 
2927  ! Store the file index a bit easier and set the base to 1.
2928 
2929  cgnsind = fileids(solid)
2930  cgnsbase = 1
2931 
2932  ! Determine the number of zones (blocks) in the restart file
2933  ! and check if this is identical to the number in the grid file.
2934 
2935  call cg_nzones_f(cgnsind, cgnsbase, nzones, ierr)
2936  if (ierr /= all_ok) &
2937  call terminate("readRestartFile", &
2938  "Something wrong when calling cg_nzones_f")
2939 
2940  if (nzones /= cgnsndom) &
2941  call terminate("readRestartFile", &
2942  "Number of blocks in grid file and restart &
2943  &file differ")
2944 
2945  ! Create a sorted version of the zone names of the restart file
2946  ! and store its corresponding zone numbers in zoneNumbers.
2947 
2949 
2950  ! Loop over the number of blocks stored on this processor.
2951 
2952  domains: do nn = 1, ndom
2953 
2954  ! Set the pointers for this block. Make sure that the
2955  ! correct data is set.
2956 
2957  ii = min(solid, ntimeintervalsspectral)
2958  call setpointers(nn, 1_inttype, ii)
2959 
2960  ! Store the zone name of the original grid a bit easier.
2961 
2962  cgnsname = cgnsdoms(nbkglobal)%zoneName
2963 
2964  ! Search in the sorted zone names of the restart file for
2965  ! cgnsName. The name must be found; otherwise the restart
2966  ! is pointless. If found, the zone number is set accordingly.
2967 
2968  jj = bsearchstrings(cgnsname, zonenames)
2969  if (jj == 0) then
2970  write (errormessage, *) "Zone name ", trim(cgnsname), &
2971  " not found in restart file ", &
2972  trim(solfiles(solid))
2973  call terminate("readRestartFile", errormessage)
2974  else
2975  jj = zonenumbers(jj)
2976  end if
2977 
2978  cgnszone = jj
2979 
2980  ! Determine the dimensions of the zone and check if these are
2981  ! identical to the dimensions of the block in the grid file.
2982 
2983  call cg_zone_read_f(cgnsind, cgnsbase, cgnszone, &
2984  cgnsname, sizes, ierr)
2985  if (ierr /= all_ok) &
2986  call terminate("readRestartFile", &
2987  "Something wrong when calling &
2988  &cg_zone_read_f")
2989 
2990  if (cgnsdoms(nbkglobal)%il /= sizes(1) .or. &
2991  cgnsdoms(nbkglobal)%jl /= sizes(2) .or. &
2992  cgnsdoms(nbkglobal)%kl /= sizes(3)) &
2993  call terminate("readRestartFile", &
2994  "Corresponding zones in restart file and &
2995  &grid file have different dimensions")
2996 
2997  ! Determine the number of flow solutions in this zone and
2998  ! check if there is a solution stored.
2999 
3000  call cg_nsols_f(cgnsind, cgnsbase, cgnszone, nsols, ierr)
3001  if (ierr /= all_ok) &
3002  call terminate("readRestartFile", &
3003  "Something wrong when calling cg_nsols_f")
3004 
3005  if (nsols == 0) &
3006  call terminate("readRestartFile", &
3007  "No solution present in restart file")
3008 
3009  ! Check for multiple solutions. A distinction is needed for
3010  ! overset cases because there will be an extra solution node
3011  ! for the nodal iblanks.
3012 
3013  if ((nsols > 1) .or. nsols > 2) &
3014  call terminate("readRestartFile", &
3015  "Multiple solutions present in restart file")
3016 
3017  ! Determine the location of the solution variables. A loop is
3018  ! done over the solution nodes which is either 1 or 2. In the
3019  ! latter case, pick the node not named "Nodal Blanks".
3020 
3021  do cgnssol = 1, nsols
3022  call cg_sol_info_f(cgnsind, cgnsbase, cgnszone, cgnssol, &
3023  cgnsname, location, ierr)
3024  if (ierr /= all_ok) &
3025  call terminate("readRestartFile", &
3026  "Something wrong when calling &
3027  &cg_sol_info_f")
3028 
3029  if (trim(cgnsname) /= "Nodal Blanks") exit
3030  end do
3031 
3032  ! Determine the rind info.
3033 
3034  call cg_goto_f(cgnsind, cgnsbase, ierr, "Zone_t", &
3035  cgnszone, "FlowSolution_t", cgnssol, "end")
3036  if (ierr /= all_ok) &
3037  call terminate("readRestartFile", &
3038  "Something wrong when calling cg_goto_f")
3039 
3040  call cg_rind_read_f(rindsizes, ierr)
3041  if (ierr /= all_ok) &
3042  call terminate("readRestartFile", &
3043  "Something wrong when calling &
3044  &cg_rind_read_f")
3045 
3046  ! Check if halo's are present. If not, set halosRead to .false.
3047  ! This only needs to be done if this is not an older state for
3048  ! an unsteady computation.
3049 
3050  if (solid == 1 .or. equationmode == timespectral) then
3051  if (rindsizes(1) == 0 .or. rindsizes(2) == 0 .or. rindsizes(3) == 0 .or. &
3052  rindsizes(4) == 0 .or. rindsizes(5) == 0 .or. rindsizes(6) == 0) &
3053  halosread = .false.
3054  end if
3055 
3056  ! Initialize the number of halo cells to read to 0.
3057 
3058  nhimin = 0; nhjmin = 0; nhkmin = 0
3059  nhimax = 0; nhjmax = 0; nhkmax = 0
3060 
3061  ! Determine the range which must be read. A few things must be
3062  ! taken into account: - in iBegor, iEndor, etc. The nodal
3063  ! range is stored. As in CGNS the cell
3064  ! range start at 1, 1 must be subtracted
3065  ! from the upper bound.
3066  ! - the rind info must be taken into
3067  ! account, because only the upper bound
3068  ! is changed in cgns; the lower bound
3069  ! remains 1.
3070  ! - in case the solution is stored in the
3071  ! vertices one extra variable in each
3072  ! direction is read. An averaging will
3073  ! take place to obtain cell centered
3074  ! values.
3075  ! Also when vertex data is present, set halosRead to .false.,
3076  ! because it is not possible to determine the halos.
3077 
3078  if (location == cellcenter) then
3079 
3080  ! Correct the number of halo cells to be read.
3081  ! Only if this is not an older state in time for an
3082  ! unsteady computation.
3083 
3084  if (solid == 1 .or. equationmode == timespectral) then
3085  if (rindsizes(1) > 0) nhimin = 1; if (rindsizes(2) > 0) nhimax = 1
3086  if (rindsizes(3) > 0) nhjmin = 1; if (rindsizes(4) > 0) nhjmax = 1
3087  if (rindsizes(5) > 0) nhkmin = 1; if (rindsizes(6) > 0) nhkmax = 1
3088  end if
3089 
3090  ! Set the cell range to be read from the CGNS file.
3091 
3092  rangemin(1) = ibegor - nhimin
3093  rangemin(2) = jbegor - nhjmin
3094  rangemin(3) = kbegor - nhkmin
3095 
3096  rangemax(1) = rangemin(1) + nx - 1 + nhimin + nhimax
3097  rangemax(2) = rangemin(2) + ny - 1 + nhjmin + nhjmax
3098  rangemax(3) = rangemin(3) + nz - 1 + nhkmin + nhkmax
3099 
3100  else if (location == vertex) then
3101 
3102  ! Set the nodal range such that enough info is present
3103  ! to average the nodal data to create the cell centered
3104  ! data in the owned cells. No halo cells will be
3105  ! initialized.
3106 
3107  halosread = .false.
3108 
3109  rangemin(1) = ibegor
3110  rangemin(2) = jbegor
3111  rangemin(3) = kbegor
3112 
3113  rangemax(1) = rangemin(1) + nx
3114  rangemax(2) = rangemin(2) + ny
3115  rangemax(3) = rangemin(3) + nz
3116  else
3117  call terminate("readRestartFile", &
3118  "Only CellCenter or Vertex data allowed in &
3119  &restart file")
3120  end if
3121 
3122  ! Allocate the memory for buffer, needed to store the variable
3123  ! to be read, and bufferVertex in case the solution is stored
3124  ! in the vertices.
3125 
3126  allocate (buffer(2 - nhimin:il + nhimax, &
3127  2 - nhjmin:jl + nhjmax, &
3128  2 - nhkmin:kl + nhkmax), stat=ierr)
3129  if (ierr /= 0) &
3130  call terminate("readRestartFile", &
3131  "Memory allocation failure for buffer")
3132 
3133  if (location == vertex) then
3134  allocate (buffervertex(1:il, 1:jl, 1:kl), stat=ierr)
3135  if (ierr /= 0) &
3136  call terminate("readRestartFile", &
3137  "Memory allocation failure for bufferVertex")
3138  end if
3139 
3140  ! Create a sorted version of the variable names and store the
3141  ! corresponding type in varTypes.
3142 
3143  call getsortedvarnumbers
3144 
3145  ! Read the density and the turbulence variables.
3146 
3147  call readdensity(ntypemismatch)
3148  call readturbvar(ntypemismatch)
3149 
3150  ! Read the other variables, depending on the situation.
3151 
3152  testprim: if (solid == 1 .or. equationmode == timespectral) then
3153 
3154  ! Either the first solution or time spectral mode. Read
3155  ! the primitive variables from the restart file.
3156 
3157  call readxvelocity(ntypemismatch)
3158  call readyvelocity(ntypemismatch)
3159  call readzvelocity(ntypemismatch)
3160  call readpressure(ntypemismatch)
3161 
3162  else testprim
3163 
3164  ! Old solution in unsteady mode. Read the conservative
3165  ! variables.
3166 
3167  call readxmomentum(ntypemismatch)
3168  call readymomentum(ntypemismatch)
3169  call readzmomentum(ntypemismatch)
3170  call readenergy(ntypemismatch)
3171 
3172  end if testprim
3173 
3174  ! Release the memory of buffer, varNames and varTypes.
3175 
3176  deallocate (buffer, varnames, vartypes, stat=ierr)
3177  if (ierr /= 0) &
3178  call terminate("readRestartFile", &
3179  "Deallocation error for buffer, varNames &
3180  &and varTypes.")
3181 
3182  ! In case bufferVertex is allocated, release it.
3183 
3184  if (location == vertex) then
3185  deallocate (buffervertex, stat=ierr)
3186  if (ierr /= 0) &
3187  call terminate("readRestartFile", &
3188  "Deallocation error for bufferVertex")
3189  end if
3190 
3191  end do domains
3192 
3193  ! Release the memory of zoneNames and zoneNumbers.
3194 
3195  deallocate (zonenames, zonenumbers, stat=ierr)
3196  if (ierr /= 0) &
3197  call terminate("readRestartFile", &
3198  "Deallocation failure for zoneNames &
3199  &and zoneNumbers.")
3200 
3201  ! Close the cgns solution file.
3202 
3203  call cg_close_f(cgnsind, ierr)
3204  if (ierr /= all_ok) &
3205  call terminate("readRestartFile", &
3206  "Something wrong when calling cg_close_f")
3207 
3208  end do solloop
3209 
3210  ! Determine the global sum of nTypeMismatch; the result only
3211  ! needs to be known on processor 0. Use ii as the global buffer
3212  ! to store the result. If a type mismatch occured,
3213  ! print a warning.
3214 
3215  call mpi_reduce(ntypemismatch, ii, 1, adflow_integer, &
3216  mpi_sum, 0, adflow_comm_world, ierr)
3217  if (myid == 0 .and. ii > 0) then
3218 
3219  write (integerstring, "(i6)") ii
3220  integerstring = adjustl(integerstring)
3221 
3222  print "(a)", "#"
3223  print "(a)", "# Warning"
3224  print strings, "# ", trim(integerstring), " type mismatches occured when reading the solution of the blocks"
3225  print "(a)", "#"
3226  end if
3227 
3228  end subroutine readrestartfile
3229 
3231  !
3232  ! getSortedZoneNumbers reads the names of the zones of the
3233  ! cgns file given by cgnsInd and cgnsBase. Afterwards the
3234  ! zonenames are sorted in increasing order, such that a binary
3235  ! search algorithm can be employed. The original zone numbers
3236  ! are stored in zoneNumbers.
3237  ! If the zone contains a link to a zone containing the
3238  ! coordinates the name of the linked zone is taken.
3239  !
3240 
3241  use constants
3242  use cgnsgrid, only: cgnsndom
3243  use su_cgns
3245  use sorting, only: qsortstrings, bsearchstrings
3246  use utils, only: terminate
3247  use commonformats, only: strings
3248  implicit none
3249  !
3250  ! Local variables.
3251  !
3252  integer :: ierr
3253  integer :: zone, zonetype, ncoords, pathlength
3254  integer :: pos
3255 
3256  integer(kind=cgsize_t), dimension(9) :: sizesblock
3257 
3258  integer(kind=intType) :: nn, ii
3259 
3260  character(len=maxStringLen) :: errormessage, linkpath
3261  character(len=maxCGNSNameLen), dimension(cgnsNdom) :: tmpnames
3262 
3263  logical :: namefound
3264 
3265  character(len=7) :: int1string, int2string
3266 
3267  ! Allocate the memory for zoneNames and zoneNumbers.
3268 
3269  allocate (zonenames(cgnsndom), zonenumbers(cgnsndom), stat=ierr)
3270  if (ierr /= 0) &
3271  call terminate("getSortedZoneNumbers", &
3272  "Memory allocation failure for zoneNames &
3273  &and zoneNumbers")
3274 
3275  ! Loop over the number of zones in the file.
3276 
3277  cgnsdomains: do nn = 1, cgnsndom
3278 
3279  ! Initialize nameFound to .false.
3280 
3281  namefound = .false.
3282 
3283  ! Check if the zone is structured.
3284 
3285  zone = nn
3286  call cg_zone_type_f(cgnsind, cgnsbase, zone, zonetype, ierr)
3287  if (ierr /= all_ok) &
3288  call terminate("getSortedZoneNumbers", &
3289  "Something wrong when calling cg_zone_type_f")
3290 
3291  if (zonetype /= structured) then
3292 
3293  write (int1string, "(i7)") cgnsbase
3294  int1string = adjustl(int1string)
3295  write (int2string, "(i7)") zone
3296  int2string = adjustl(int2string)
3297 
3298  write (errormessage, strings) "Base ", trim(int1string), ": Zone ", trim(int2string), &
3299  " of the cgns restart file is not structured"
3300  call terminate("getSortedZoneNumbers", errormessage)
3301 
3302  end if
3303 
3304  ! Determine the number of grid coordinates of this zone.
3305 
3306  call cg_ncoords_f(cgnsind, cgnsbase, zone, ncoords, ierr)
3307  if (ierr /= all_ok) &
3308  call terminate("getSortedZoneNumbers", &
3309  "Something wrong when calling cg_ncoords_f")
3310 
3311  ! If ncoords == 3, there are coordinates present. Then check if
3312  ! it is a link. If so, take the zone name of the link.
3313 
3314  if (ncoords == 3) then
3315 
3316  ! Go to the coordinates node.
3317 
3318  call cg_goto_f(cgnsind, cgnsbase, ierr, "Zone_t", zone, &
3319  "GridCoordinates_t", 1, "end")
3320  if (ierr /= all_ok) &
3321  call terminate("getSortedZoneNumbers", &
3322  "Something wrong when calling cg_goto_f")
3323 
3324  ! Check if this node is a link.
3325 
3326  call cg_is_link_f(pathlength, ierr)
3327  if (ierr /= all_ok) &
3328  call terminate("getSortedZoneNumbers", &
3329  "Something wrong when calling cg_is_link_f")
3330 
3331  if (pathlength > 0) then
3332 
3333  ! Determine the name of the linkPath.
3334 
3335  call cg_link_read_f(errormessage, linkpath, ierr)
3336  if (ierr /= all_ok) &
3337  call terminate("getSortedZoneNumbers", &
3338  "Something wrong when calling &
3339  &cg_link_read_f")
3340 
3341  ! Find the zone name.
3342  ! Find, starting from the back, the forward slash.
3343 
3344  pos = index(linkpath, "/", .true.)
3345  if (pos > 0) then
3346  linkpath = linkpath(:pos - 1)
3347 
3348  ! Find the next forward slash from the back and
3349  ! remove the leading part from the path name.
3350 
3351  pos = index(linkpath, "/", .true.)
3352  if (pos > 0) linkpath = linkpath(pos + 1:)
3353  end if
3354 
3355  ! Create the zone name and set nameFound to .true..
3356 
3357  linkpath = adjustl(linkpath)
3358  zonenames(nn) = trim(linkpath)
3359  namefound = .true.
3360 
3361  end if
3362  end if
3363 
3364  ! If no name was yet found set it to the name of the
3365  ! current zone.
3366 
3367  if (.not. namefound) then
3368  call cg_zone_read_f(cgnsind, cgnsbase, zone, &
3369  zonenames(nn), sizesblock, ierr)
3370  if (ierr /= all_ok) &
3371  call terminate("getSortedZoneNumbers", &
3372  "Something wrong when calling &
3373  &cg_zone_read_f")
3374  end if
3375 
3376  end do cgnsdomains
3377 
3378  ! Set tmpNames to zoneNames and sort the latter
3379  ! in increasing order.
3380 
3381  do nn = 1, cgnsndom
3382  tmpnames(nn) = zonenames(nn)
3383  end do
3384 
3385  ! Sort zoneNames in increasing order.
3386 
3388 
3389  ! Initialize zoneNumbers to -1. This serves as a check during
3390  ! the search.
3391 
3392  do nn = 1, cgnsndom
3393  zonenumbers(nn) = -1
3394  end do
3395 
3396  ! Find the original zone numbers for the sorted zone names.
3397 
3398  do nn = 1, cgnsndom
3399  ii = bsearchstrings(tmpnames(nn), zonenames)
3400 
3401  ! Check if the zone number is not already taken. If this is the
3402  ! case, this means that two identical zone names are present.
3403 
3404  if (zonenumbers(ii) /= -1) &
3405  call terminate("getSortedZoneNumbers", &
3406  "Error occurs only when two identical zone &
3407  &names are present")
3408 
3409  ! And set the zone number.
3410 
3411  zonenumbers(ii) = nn
3412  end do
3413 
3414  end subroutine getsortedzonenumbers
3415 
3417  !
3418  ! getSortedVarNumbers reads the names of variables stored in
3419  ! the given solution node of the cgns file, indicated by
3420  ! cgnsInd, cgnsBase and cgnsZone. Afterwards the variable
3421  ! names are sorted in increasing order, such that they can be
3422  ! used in a binary search. Their original variable number and
3423  ! type is stored.
3424  !
3425  use constants
3426  use su_cgns
3429  use sorting, only: qsortstrings, bsearchstrings
3430  use utils, only: terminate
3431  implicit none
3432  !
3433  ! Local variables.
3434  !
3435  integer :: i, ierr
3436  integer, dimension(:), allocatable :: tmptypes
3437 
3438  integer(kind=intType) :: nn, ii
3439 
3440  integer(kind=intType), dimension(:), allocatable :: varnumbers
3441 
3442  character(len=maxCGNSNameLen), allocatable, dimension(:) :: &
3443  tmpnames
3444 
3445  ! Determine the number of solution variables stored.
3446 
3447  call cg_nfields_f(cgnsind, cgnsbase, cgnszone, cgnssol, &
3448  nvar, ierr)
3449  if (ierr /= all_ok) &
3450  call terminate("getSortedVarNumbers", &
3451  "Something wrong when calling cg_nfield_f")
3452 
3453  ! Allocate the memory for varnames, vartypes and varnumber
3454 
3455  allocate (varnames(nvar), vartypes(nvar), varnumbers(nvar), &
3456  stat=ierr)
3457  if (ierr /= 0) &
3458  call terminate("getSortedVarNumbers", &
3459  "Memory allocation failure for varNames, etc.")
3460 
3461  ! Loop over the number of variables and store their names and
3462  ! types.
3463 
3464  do i = 1, nvar
3465  call cg_field_info_f(cgnsind, cgnsbase, cgnszone, cgnssol, &
3466  i, vartypes(i), varnames(i), ierr)
3467  if (ierr /= 0) &
3468  call terminate("getSortedVarNumbers", &
3469  "Something wrong when calling cg_field_info_f")
3470  end do
3471 
3472  ! Allocate the memory for tmpTypes and tmpNames and initialize
3473  ! their values.
3474 
3475  allocate (tmptypes(nvar), tmpnames(nvar), stat=ierr)
3476  if (ierr /= 0) &
3477  call terminate("getSortedVarNumbers", &
3478  "Memory allocation failure for tmp variables")
3479 
3480  do i = 1, nvar
3481  tmptypes(i) = vartypes(i)
3482  tmpnames(i) = varnames(i)
3483  end do
3484 
3485  ! Sort varNames in increasing order.
3486 
3487  nn = nvar
3488  call qsortstrings(varnames, nn)
3489 
3490  ! Initialize varNumbers to -1. This serves as a check during
3491  ! the search.
3492 
3493  do i = 1, nvar
3494  varnumbers(i) = -1
3495  end do
3496 
3497  ! Find the original types and numbers for the just sorted
3498  ! variable names.
3499 
3500  do i = 1, nvar
3501  ii = bsearchstrings(tmpnames(i), varnames)
3502 
3503  ! Check if the variable number is not already taken. If this is
3504  ! the case, this means that two identical var names are present.
3505 
3506  if (varnumbers(ii) /= -1) &
3507  call terminate("getSortedVarNumbers", &
3508  "Error occurs only when two identical &
3509  &variable names are present")
3510 
3511  ! And set the variable number and type.
3512 
3513  varnumbers(ii) = i
3514  vartypes(ii) = tmptypes(i)
3515  end do
3516 
3517  ! Release the memory of varNumbers, tmpNames and tmpTypes.
3518 
3519  deallocate (varnumbers, tmptypes, tmpnames, stat=ierr)
3520  if (ierr /= 0) &
3521  call terminate("getSortedVarNumbers", &
3522  "Deallocation error for tmp variables")
3523 
3524  end subroutine getsortedvarnumbers
3525 #endif
3526 end module initializeflow
Definition: BCData.F90:1
subroutine setbcdatacoarsegrid
Definition: BCData.F90:3091
subroutine setbcdatafinegrid(initializationPart)
Definition: BCData.F90:2623
subroutine applyallbc_block(secondHalo)
Definition: BCRoutines.F90:58
subroutine applyallbc(secondHalo)
Definition: BCRoutines.F90:16
Definition: block.F90:1
integer(kind=inttype) ndom
Definition: block.F90:761
type(blocktype), dimension(:, :, :), allocatable, target flowdoms
Definition: block.F90:771
integer(kind=inttype) jl
integer(kind=inttype) kbegor
integer(kind=inttype) nx
real(kind=realtype), dimension(:, :, :), pointer p
integer(kind=inttype) ny
integer(kind=inttype) ie
real(kind=realtype), dimension(:, :, :, :), pointer w
real(kind=realtype), dimension(:, :, :), pointer d2wall
integer(kind=inttype) nbkglobal
integer(kind=inttype) jb
integer(kind=inttype) kb
real(kind=realtype), dimension(:, :, :), pointer rlv
integer(kind=inttype) ibegor
integer(kind=inttype) nbocos
integer(kind=inttype) sectionid
integer(kind=inttype) jbegor
real(kind=realtype), dimension(:, :, :), pointer rev
integer(kind=inttype), dimension(:), pointer bctype
real(kind=realtype), dimension(:, :, :, :), pointer dw
integer(kind=inttype) ib
real(kind=realtype), dimension(:, :, :, :), pointer fw
integer(kind=inttype) nz
integer(kind=inttype) je
integer(kind=inttype) ke
integer(kind=inttype) kl
integer(kind=inttype) il
type(cgnsblockinfotype), dimension(:), allocatable cgnsdoms
Definition: cgnsGrid.F90:495
integer(kind=inttype) cgnsndom
Definition: cgnsGrid.F90:491
character(len=maxstringlen) strings
character(len=maxstringlen) stringint1
integer adflow_comm_world
integer(kind=inttype), parameter spalartallmarasedwards
Definition: constants.F90:128
integer(kind=inttype), parameter spalartallmaras
Definition: constants.F90:128
real(kind=realtype), parameter zero
Definition: constants.F90:71
integer(kind=inttype), parameter coupled
Definition: constants.F90:207
integer, parameter itu3
Definition: constants.F90:43
real(kind=realtype), parameter third
Definition: constants.F90:81
real(kind=realtype), parameter pi
Definition: constants.F90:22
integer, parameter itu2
Definition: constants.F90:42
real(kind=realtype), parameter eps
Definition: constants.F90:23
integer(kind=inttype), parameter ktau
Definition: constants.F90:128
integer(kind=inttype), parameter md
Definition: constants.F90:189
integer, parameter itu1
Definition: constants.F90:40
integer, parameter irho
Definition: constants.F90:34
integer(kind=inttype), parameter komegawilcox
Definition: constants.F90:128
integer(kind=inttype), parameter timespectral
Definition: constants.F90:115
integer(kind=inttype), parameter komegamodified
Definition: constants.F90:128
integer, parameter ivx
Definition: constants.F90:35
integer(kind=inttype), parameter unsteady
Definition: constants.F90:115
integer, parameter irhoe
Definition: constants.F90:38
integer(kind=inttype), parameter bdf
Definition: constants.F90:189
integer(kind=inttype), parameter totalconditions
Definition: constants.F90:235
real(kind=realtype), parameter one
Definition: constants.F90:72
integer, parameter maxstringlen
Definition: constants.F90:16
integer(kind=inttype), parameter steady
Definition: constants.F90:115
integer, parameter ivz
Definition: constants.F90:37
real(kind=realtype), parameter two
Definition: constants.F90:73
integer, parameter itu4
Definition: constants.F90:44
integer(kind=inttype), parameter explicitrk
Definition: constants.F90:189
integer(kind=inttype), parameter massflow
Definition: constants.F90:235
integer(kind=inttype), parameter mentersst
Definition: constants.F90:128
integer(kind=inttype), parameter subsonicinflow
Definition: constants.F90:265
integer, parameter ivy
Definition: constants.F90:36
integer(kind=inttype), parameter ransequations
Definition: constants.F90:110
integer(kind=inttype), parameter internalflow
Definition: constants.F90:120
integer(kind=inttype), parameter v2f
Definition: constants.F90:128
subroutine computelamviscosity(includeHalos)
Definition: flowUtils.F90:1202
subroutine adjustinflowangle()
Definition: flowUtils.F90:1326
subroutine computeetotblock(iStart, iEnd, jStart, jEnd, kStart, kEnd, correctForK)
Definition: flowUtils.F90:553
subroutine computegamma(T, gamma, mm)
Definition: flowUtils.F90:56
subroutine etot(rho, u, v, w, p, k, etotal, correctForK)
Definition: flowUtils.F90:675
real(kind=realtype) rhoinfdim
real(kind=realtype) muinfdim
real(kind=realtype) gammainf
real(kind=realtype) href
integer(kind=inttype) nt1
real(kind=realtype) uinf
real(kind=realtype) pinfdim
real(kind=realtype) pinfcorr
real(kind=realtype) pinf
real(kind=realtype) rgas
real(kind=realtype) uref
real(kind=realtype) muref
real(kind=realtype) tinfdim
real(kind=realtype) rhoref
integer(kind=inttype) nwf
real(kind=realtype), dimension(:), allocatable winf
real(kind=realtype) muinf
real(kind=realtype) tref
integer(kind=inttype) nw
real(kind=realtype) rhoinf
real(kind=realtype) pref
real(kind=realtype) timeref
integer(kind=inttype) nt2
subroutine whalo2(level, start, end, commPressure, commGamma, commViscous)
subroutine allocmemflovarpart2(sps, level)
subroutine infchangecorrection(oldWinf, correctionTol, correctionType)
subroutine getsortedvarnumbers
subroutine updatebcdataalllevels
subroutine setpressureandcomputeenergy(halosRead)
subroutine allocmemflovarpart1(sps, level)
subroutine determinesolfilenames
subroutine velmagnanddirectionsubface(vmag, dir, BCData, mm)
subroutine setrestartfiles(fileName, i)
subroutine initdepvarandhalos(halosRead)
subroutine setsolfilenames
subroutine interpolatespectralsolution
subroutine copyspectralsolution
subroutine setuniformflow
subroutine initflow
subroutine checksolfilenames
subroutine getsortedzonenumbers
subroutine timespectralcoef(coefSpectral, matrixCoefSpectral, diagMatCoefSpectral)
subroutine initializehalos(halosRead)
subroutine timespectralmatrices
subroutine referencestate
subroutine releaseextramembcs
subroutine initflowrestart
subroutine allocrestartfiles(nFiles)
subroutine readrestartfile()
subroutine initflowfield
character(len=maxstringlen), dimension(:), allocatable restartfiles
Definition: inputParam.F90:161
integer(kind=inttype) ncycles
Definition: inputParam.F90:262
integer(kind=inttype) mgstartlevel
Definition: inputParam.F90:270
integer(kind=inttype) turbtreatment
Definition: inputParam.F90:269
integer(kind=inttype) nsgstartup
Definition: inputParam.F90:264
real(kind=realtype) ssuthdim
Definition: inputParam.F90:606
real(kind=realtype) turbintensityinf
Definition: inputParam.F90:599
real(kind=realtype) gammaconstant
Definition: inputParam.F90:597
real(kind=realtype) eddyvisinfratio
Definition: inputParam.F90:599
integer(kind=inttype) equations
Definition: inputParam.F90:585
integer(kind=inttype) equationmode
Definition: inputParam.F90:585
real(kind=realtype) tsuthdim
Definition: inputParam.F90:606
integer(kind=inttype) turbmodel
Definition: inputParam.F90:586
real(kind=realtype) machcoef
Definition: inputParam.F90:595
integer(kind=inttype) flowtype
Definition: inputParam.F90:585
real(kind=realtype) mach
Definition: inputParam.F90:595
real(kind=realtype) rgasdim
Definition: inputParam.F90:597
real(kind=realtype), dimension(3) veldirfreestream
Definition: inputParam.F90:601
real(kind=realtype) musuthdim
Definition: inputParam.F90:606
real(kind=realtype), dimension(:, :, :), allocatable dvector
Definition: inputParam.F90:664
real(kind=realtype), dimension(:, :, :), allocatable rotmatrixspectral
Definition: inputParam.F90:699
integer(kind=inttype) ntimeintervalsspectral
Definition: inputParam.F90:649
real(kind=realtype), dimension(:, :, :), allocatable dscalar
Definition: inputParam.F90:663
logical useale
Definition: inputParam.F90:758
integer(kind=inttype) nrkstagesunsteady
Definition: inputParam.F90:746
integer(kind=inttype) timeintegrationscheme
Definition: inputParam.F90:723
integer(kind=inttype) ntimestepsfine
Definition: inputParam.F90:735
type(iotype), dimension(:, :), allocatable iovar
Definition: IOModule.f90:49
integer(kind=inttype) noldlevels
Definition: iteration.f90:79
integer(kind=inttype) currentlevel
Definition: iteration.f90:18
integer(kind=inttype) groundlevel
Definition: iteration.f90:18
integer(kind=inttype) nalesteps
Definition: iteration.f90:107
logical, dimension(:), allocatable oldsolwritten
Definition: iteration.f90:123
integer(kind=inttype) noldsolavail
Definition: iteration.f90:79
real(kind=realtype) timeunsteadyrestart
Definition: monitor.f90:98
integer(kind=inttype) ntimestepsrestart
Definition: monitor.f90:97
logical oversetpresent
Definition: overset.F90:373
subroutine residual
Definition: residuals.F90:1029
integer(kind=inttype) nsections
Definition: section.f90:44
type(sectiontype), dimension(:), allocatable sections
Definition: section.f90:46
subroutine gridvelocitiescoarselevels(sps)
subroutine slipvelocitiesfinelevel(useOldCoor, t, sps)
subroutine slipvelocitiescoarselevels(sps)
subroutine normalvelocitiesalllevels(sps)
subroutine gridvelocitiesfinelevel(useOldCoor, t, sps)
integer(kind=inttype) function bsearchstrings(key, base)
Definition: sorting.F90:812
subroutine qsortstrings(arr, nn)
Definition: sorting.F90:525
subroutine bcturbtreatment
subroutine applyallturbbc(secondHalo)
subroutine applyallturbbcthisblock(secondHalo)
real(kind=realtype) function sanuknowneddyratio(eddyRatio, nuLam)
Definition: turbUtils.F90:334
subroutine computeeddyviscosity(includeHalos)
Definition: turbUtils.F90:582
Definition: utils.F90:1
subroutine alloctimearrays(nTimeTot)
Definition: utils.F90:5886
subroutine rotmatrixrigidbody(tNew, tOld, rotationMatrix, rotationPoint)
Definition: utils.F90:614
real(kind=realtype) function mynorm2(x)
Definition: utils.F90:1697
subroutine spectralinterpolcoef(nsps, t, alpScal, alpMat)
Definition: utils.F90:3668
subroutine setpointers(nn, mm, ll)
Definition: utils.F90:3237
subroutine allocconvarrays(nIterTot)
Definition: utils.F90:5830
subroutine getrotmatrix(vec1, vec2, rotMat)
Definition: utils.F90:6592
subroutine terminate(routineName, errorMessage)
Definition: utils.F90:502
real(kind=cgnsrealtype), dimension(:, :, :), allocatable buffer
subroutine readymomentum(nTypeMismatch)
integer(kind=inttype) nsolsread
subroutine readyvelocity(nTypeMismatch)
integer(kind=inttype) solid
subroutine readturbvar(nTypeMismatch)
integer, dimension(:), allocatable vartypes
subroutine readzvelocity(nTypeMismatch)
integer(kind=inttype), dimension(:), allocatable zonenumbers
real(kind=cgnsrealtype), dimension(:, :, :), allocatable buffervertex
subroutine readdensity(nTypeMismatch)
subroutine readzmomentum(nTypeMismatch)
character(len=maxstringlen), dimension(:), allocatable solfiles
subroutine readxvelocity(nTypeMismatch)
subroutine readpressure(nTypeMismatch)
subroutine readenergy(nTypeMismatch)
character(len=maxcgnsnamelen), dimension(:), allocatable zonenames
subroutine readxmomentum(nTypeMismatch)
integer(kind=cgsize_t), dimension(3) rangemax
character(len=maxcgnsnamelen), dimension(:), allocatable varnames
integer(kind=cgsize_t), dimension(3) rangemin
subroutine scalefactors(fileIDs)