ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
block.F90
Go to the documentation of this file.
1 module block
2  !
3  ! This module contains the definition of the derived data type
4  ! for block, which is the basic building block for this code.
5  ! Apart from the derived data type for block, this module also
6  ! contains the actual array for storing the blocks and the
7  ! number of blocks stored on this processor.
8  !
9  use constants, only: realtype, inttype, portype, maxcgnsnamelen, &
11  implicit none
12  save
13 
14  !
15  ! The definition of the derived data type visc_subface_type,
16  ! which stores the viscous stress tensor and heat flux vector.
17  ! In this way it is avoided that these quantities must be
18  ! recomputed for the viscous forces and postprocessing. This
19  ! saves both time and a considerable amount of code.
20  !
22 
23  ! tau(:,:,6): The 6 components of the viscous stress tensor.
24  ! The first 2 dimensions of these arrays are equal
25  ! to the dimenions of the cell subface without any
26  ! halo cell. Consequently the starting index is
27  ! arbitrary, such that no offset computation is
28  ! needed when the arrays are accessed.
29  ! q(:,:,3): Same story for the heat flux vector.
30  ! uTau(:,:): And for the friction velocity.
31 
32  real(kind=realtype), dimension(:, :, :), pointer :: tau, q
33  real(kind=realtype), dimension(:, :), pointer :: utau
34 
35  end type viscsubfacetype
36 
37  type rptr
38  real(kind=realtype), dimension(:, :, :), pointer :: var
39  end type rptr
40 
41  type iptr
42  integer(kind=intType), dimension(:, :, :), pointer :: var
43  end type iptr
44  !
45  ! The definition of the derived data type BCDataType, which
46  ! stores the prescribed data of boundary faces as well as unit
47  ! normals. For all the arrays the first two dimensions equal the
48  ! dimensions of the subface, possibly extended with halo cells.
49  ! Consequently the starting index is arbitrary, such that no
50  ! offset computation is needed when the array is accessed.
51  !
53 
54  ! inBeg, inEnd: Node range in the first direction of the subface
55  ! jnBeg, jnEnd: Idem in the second direction.
56  ! icBeg, icEnd: Cell range in the first direction of the subface
57  ! jcBeg, jcEnd: Idem in the second direction.
58 
59  integer(kind=intType) :: inbeg, inend, jnbeg, jnend
60  integer(kind=intType) :: icbeg, icend, jcbeg, jcend
61 
62  ! norm(:,:,3): The unit normal; it points out of the domain.
63  ! rface(:,:): Velocity of the face in the direction of the
64  ! outward pointing normal. only allocated for
65  ! the boundary conditions that need this info.
66 
67  real(kind=realtype), dimension(:, :, :), pointer :: norm
68  real(kind=realtype), dimension(:, :), pointer :: rface
69  real(kind=realtype), dimension(:, :, :), pointer :: f, fv, fp
70  real(kind=realtype), dimension(:, :, :), pointer :: t, tv, tp
71  real(kind=realtype), dimension(:, :), pointer :: area
72  real(kind=realtype), dimension(:, :), pointer :: cptarget
73  integer(kind=realType), dimension(:, :), pointer :: surfindex
74 
75  ! Generic pointers for performing a globalized reduction.
76  real(kind=realtype), dimension(:, :), pointer :: nodeval
77  real(kind=realtype), dimension(:, :), pointer :: cellval
78 
79  ! symNorm is the normal for (symmertry) boundary conditions.
80  ! symNormSet is set to false until symNorm is computed at the
81  ! beginning of a simulation. symNorm then remains constant for
82  ! the remainder of the simulation. This is ok, since if the
83  ! normal of the symmetry plane is changing, your results are
84  ! invalid anyway. These values are only used on symmetry
85  ! plane. They are undefined for other BC's.
86  real(kind=realtype), dimension(3) :: symnorm
87 
88  logical :: symnormset
89 
90  ! subsonicInletTreatment: which boundary condition treatment
91  ! to use for subsonic inlets; either
92  ! totalConditions or massFlow.
93 
94  integer(kind=intType) :: subsonicinlettreatment
95 
96  ! uSlip(:,:,3): the 3 components of the velocity vector on
97  ! a viscous wall.
98  ! TNS_Wall(:,:): Wall temperature for isothermal walls.
99 
100  real(kind=realtype), dimension(:, :, :), pointer :: uslip
101  real(kind=realtype), dimension(:, :), pointer :: tns_wall
102 
103  ! The name of this boundary condition and it's index
104  character(maxCGNSNameLen) :: family
105  integer(kind=intType) :: famid
106 
107  ! Added by HDN
108  ! normALE(0:nALEsteps,ie:ib,je:jb,3)
109  ! - Storage of norm for intermediate meshes.
110  ! rFaceALE(0:nALEsteps,iBeg:iEnd,jBeg:jEnd)
111  ! - Storage of rface for intermediate meshes.
112  ! uSlipALE(0:nALEsteps,iBeg:iEnd,jBeg:jEnd,3)
113  ! - Storage of uSlip for intermediate meshes.
114  ! cellHeatFlux(iBeg:iEnd,jBeg:jEnd,3)
115  ! - Surface heat flux (cell based/node based).
116  real(kind=realtype), dimension(:, :, :, :), pointer :: normale
117  real(kind=realtype), dimension(:, :, :), pointer :: rfaceale
118  real(kind=realtype), dimension(:, :, :, :), pointer :: uslipale
119  real(kind=realtype), dimension(:, :), pointer :: nodeheatflux
120  real(kind=realtype), dimension(:, :), pointer :: cellheatflux
121 
122  ! ptInlet(:,:): Total pressure at subsonic inlets.
123  ! ttInlet(:,:): Total temperature at subsonic inlets.
124  ! htInlet(:,:): Total enthalpy at subsonic inlets.
125  ! flowXDirInlet(:,:): X-direction of the flow for subsonic
126  ! inlets.
127  ! flowYDirInlet(:,:): Idem in y-direction.
128  ! flowZDirInlet(:,:): Idem in z-direction.
129 
130  real(kind=realtype), dimension(:, :), pointer :: ptinlet, ttinlet, htinlet
131  real(kind=realtype), dimension(:, :), pointer :: flowxdirinlet, flowydirinlet, flowzdirinlet
132 
133  ! turbInlet(:,:,nt1:nt2): Turbulence variables at inlets,
134  ! either subsonic or supersonic.
135 
136  real(kind=realtype), dimension(:, :, :), pointer :: turbinlet
137 
138  ! rho(:,:): density; used for multiple bc's.
139  ! velX(:,:): x-velocity; used for multiple bc's.
140  ! velY(:,:): y-velocity; used for multiple bc's.
141  ! velZ(:,:): z-velocity; used for multiple bc's.
142  ! ps(:,:): static pressure; used for multiple bc's.
143 
144  real(kind=realtype), dimension(:, :), pointer :: rho
145  real(kind=realtype), dimension(:, :), pointer :: velx, vely, velz
146  real(kind=realtype), dimension(:, :), pointer :: ps
147 
148  ! Surface blanking for force integration
149  integer(kind=intType), dimension(:, :), pointer :: iblank
150 
151  ! Surface deviation. This is an estimate of how much the surface
152  ! deviates from the "real" underlying surface'
153  real(kind=realtype), dimension(:, :), pointer :: delta
154  real(kind=realtype), dimension(:, :), pointer :: deltanode
155 
156  end type bcdatatype
157 
159  real(kind=realtype), dimension(:, :, :), pointer :: weight
160  end type surfacenodeweightarray
161 
163 
164  ! Make everything in here static such that we can easily copy the
165  ! datatype and use MPI to communicate them directly. This data
166  ! type bears some resemblence to the haloList type used for the
167  ! B2B preprocessing.
168 
169  ! The quality metric for the fringe
170  real(kind=realtype) :: quality
171 
172  ! This is the information regarding where the cell came from.
173  integer(kind=intType) :: myblock, myindex
174 
175  ! This is the information about the donor that was found
176  integer(kind=intType) :: donorproc, donorblock, dindex
177  real(kind=realtype) :: donorfrac(3)
178 
179  end type fringetype
180 
181  interface operator(<=)
182  module procedure lessequalfringetype
183  end interface operator(<=)
184 
185  interface operator(<)
186  module procedure lessfringetype
187  end interface operator(<)
188 
190  integer(kind=intType) :: donorproc, donorblock, di, dj, dk, myblock
191  real(kind=realtype) :: donorfrac(3)
192  end type interppttype
193 
194  interface operator(<=)
195  module procedure lessequalinterppttype
196  end interface operator(<=)
197 
198  interface operator(<)
199  module procedure lessinterppttype
200  end interface operator(<)
201 
202  ! The definition of the derived data type block_type, which
203  ! stores dimensions, coordinates, solution, etc.
204  !
206  !
207  ! Block dimensions and orientation.
208  !
209  ! nx, ny, nz - Block integer dimensions for no halo cell based
210  ! quantities.
211  ! il, jl, kl - Block integer dimensions for no halo node based
212  ! quantities.
213  ! ie, je, ke - Block integer dimensions for single halo
214  ! cell-centered quantities.
215  ! ib, jb, kb - Block integer dimensions for double halo
216  ! cell-centered quantities.
217  ! rightHanded - Whether or not the block is a right handed.
218  ! If not right handed it is left handed of course.
219 
220  integer(kind=intType) :: nx, ny, nz
221  integer(kind=intType) :: il, jl, kl
222  integer(kind=intType) :: ie, je, ke
223  integer(kind=intType) :: ib, jb, kb
224 
225  logical :: righthanded
226  !
227  ! Block boundary conditions.
228  !
229  ! nSubface - Number of subfaces on this block.
230  ! n1to1 - Number of 1 to 1 block boundaries.
231  ! nBocos - Number of physical boundary subfaces.
232  ! nViscBocos - Number of viscous boundary subfaces.
233  ! BCType(:) - Boundary condition type for each
234  ! subface. See the module BCTypes for
235  ! the possibilities.
236  ! BCFaceID(:) - Block face location of each subface.
237  ! possible values are: iMin, iMax, jMin,
238  ! jMax, kMin, kMax. see also module
239  ! BCTypes.
240  ! cgnsSubface(:) - The subface in the corresponding cgns
241  ! block. As cgns distinguishes between
242  ! boundary and internal boundaries, the
243  ! BCType of the subface is needed to
244  ! know which one to take.
245  ! inBeg(:), inEnd(:) - Lower and upper limits for the nodes
246  ! jnBeg(:), jnEnd(:) in each of the index directions on a
247  ! knBeg(:), knEnd(:) given subface. Note that one of these
248  ! indices will not change since we will
249  ! be moving on a face.
250  ! dinBeg(:), dinEnd(:) - Lower and upper limits for the nodes
251  ! djnBeg(:), djnEnd(:) in the each of the index directions
252  ! dknBeg(:), dknEnd(:) of the donor subface for this
253  ! particular subface. Note that one of
254  ! these indices will not change since we
255  ! will be moving on a face.
256  ! icBeg(:), icEnd(:) - Lower and upper limits for the cells
257  ! jcBeg(:), jcEnd(:) in each of the index directions for
258  ! kcBeg(:), kcEnd(:) the subface. The cells indicated by
259  ! this range are halo cells (the
260  ! constant index) adjacent to the face.
261  ! a possible overlap outside the block
262  ! is stored.
263  ! neighBlock(:) - Local block number to which this
264  ! subface connects. This value is set to
265  ! zero if this subface is not connected
266  ! to another block.
267  ! neighProc(:) - Processor number where the neighbor
268  ! block is stored. This value is set to
269  ! -1 if this subface is not connected
270  ! to another block.
271  ! l1(:), l2(:), - Short hand for the transformation
272  ! l3(:) matrix between this subface and the
273  ! neighbor block. These values are set
274  ! to zero if this subface is not
275  ! connected to another block.
276  ! groupNum(:) - Group number to which this subface
277  ! belongs. If this subface does not
278  ! belong to any group, the corresponding
279  ! entry in this array is zeroed out. If
280  ! the subface belongs to a sliding mesh
281  ! interface the absolute value of
282  ! groupNum contains the number of the
283  ! sliding mesh interface. One side of
284  ! the interface gets a positive number,
285  ! the other side a negative one.
286  !
287  !
288  integer(kind=intType) :: nsubface, n1to1, nbocos, nviscbocos
289 
290  integer(kind=intType), dimension(:), pointer :: bctype
291  integer(kind=intType), dimension(:), pointer :: bcfaceid
292 
293  integer(kind=intType), dimension(:), pointer :: cgnssubface
294 
295  integer(kind=intType), dimension(:), pointer :: inbeg, inend
296  integer(kind=intType), dimension(:), pointer :: jnbeg, jnend
297  integer(kind=intType), dimension(:), pointer :: knbeg, knend
298 
299  integer(kind=intType), dimension(:), pointer :: dinbeg, dinend
300  integer(kind=intType), dimension(:), pointer :: djnbeg, djnend
301  integer(kind=intType), dimension(:), pointer :: dknbeg, dknend
302 
303  integer(kind=intType), dimension(:), pointer :: icbeg, icend
304  integer(kind=intType), dimension(:), pointer :: jcbeg, jcend
305  integer(kind=intType), dimension(:), pointer :: kcbeg, kcend
306 
307  integer(kind=intType), dimension(:), pointer :: neighblock
308  integer(kind=intType), dimension(:), pointer :: neighproc
309  integer(kind=intType), dimension(:), pointer :: l1, l2, l3
310  integer(kind=intType), dimension(:), pointer :: groupnum
311 
312  !
313  ! Overset interpolation information
314 
315  integer(kind=intType), dimension(:, :, :), pointer :: iblank
316  integer(kind=intType), dimension(:, :, :), pointer :: iblanklast
317  integer(kind=intType), dimension(:, :, :), pointer :: status
318  integer(kind=intType), dimension(:, :, :), pointer :: forcedrecv
319  type(fringetype), dimension(:), pointer :: fringes => null()
320  integer(kind=intType), dimension(:, :, :, :), pointer :: fringeptr => null()
321  integer(kind=intType), dimension(:, :, :, :), pointer :: gind => null()
322  integer(kind=intType), pointer :: ndonors
323  integer(kind=intType) :: ndonorsonownedcells
324 
325  integer(kind=intType), dimension(:, :), pointer :: orphans
326  integer(kind=intType) :: norphans
327 
328  !
329  ! Boundary data for the boundary subfaces.
330  !
331  ! BCData(nBocos): The boundary data for each of the boundary
332  ! subfaces.
333 
334  type(bcdatatype), dimension(:), pointer :: bcdata
335  !
336  ! The stress tensor and heat flux vector at viscous wall faces
337  ! as well as the face pointers to these viscous wall faces.
338  !
339  ! viscSubface(nViscBocos): Storage for the viscous stress
340  ! tensor and heat flux vector for
341  ! the viscous subfaces.
342  ! viscIMinPointer(2:jl,2:kl): Pointer to viscous subface for
343  ! the iMin block face. If the face
344  ! is not part of a viscous subface
345  ! this value is set to 0.
346  ! viscIMaxPointer(2:jl,2:kl): Idem for iMax block face.
347  ! viscJMinPointer(2:il,2:kl): Idem for jMin block face.
348  ! viscJMaxPointer(2:il,2:kl): Idem for jmax block face.
349  ! viscKMinPointer(2:il,2:jl): Idem for kMin block face.
350  ! viscKMaxPointer(2:il,2:jl): Idem for kMax block face.
351 
352  type(viscsubfacetype), dimension(:), pointer :: viscsubface
353 
354  integer(kind=intType), dimension(:, :), pointer :: visciminpointer
355  integer(kind=intType), dimension(:, :), pointer :: viscimaxpointer
356  integer(kind=intType), dimension(:, :), pointer :: viscjminpointer
357  integer(kind=intType), dimension(:, :), pointer :: viscjmaxpointer
358  integer(kind=intType), dimension(:, :), pointer :: visckminpointer
359  integer(kind=intType), dimension(:, :), pointer :: visckmaxpointer
360  !
361  ! Mesh related variables.
362  !
363  ! x(0:ie,0:je,0:ke,3) - xyz locations of grid points in block.
364  ! xOld(nOld,:,:,:,:) - Coordinates on older time levels;
365  ! only needed for unsteady problems on
366  ! deforming grids. Only allocated on
367  ! the finest grid level. The blank
368  ! dimensions are equal to the dimensions
369  ! of x.
370  ! sI(0:ie,1:je,1:ke,3) - Projected areas in the i-coordinate
371  ! direction. Normals point in the
372  ! direction of increasing i.
373  ! sJ(1:ie,0:je,1:ke,3) - Projected areas in the j-coordinate
374  ! direction. Normals point in the
375  ! direction of increasing j.
376  ! sK(1:ie,1:je,0:ke,3) - Projected areas in the k-coordinate
377  ! direction. Normals point in the
378  ! direction of increasing k.
379  ! vol(0:ib,0:jb,0:kb) - Cell volumes. The second level halo
380  ! is present for a multigrid option.
381  ! volOld(nold,2:il,..) - Volumes on older time levels; only
382  ! needed for unsteady problems on
383  ! deforming grids. Only allocated on
384  ! the finest grid level.
385  ! skew(0:ib,0:jb,0:kb) - Skewness of Volume
386  ! uv(2,2:il,2:jl,2:kl) - Parametric location on elemID for each cell.
387  ! Only used for fast wall distance calcs.
388  ! porI(1:il,2:jl,2:kl) - Porosity in the i direction.
389  ! porJ(2:il,1:jl,2:kl) - Porosity in the j direction.
390  ! porK(2:il,2:jl,1:kl) - Porosity in the k direction.
391  !
392  ! indFamilyI(:,:,:) - Index of the i-face in the arrays
393  ! to compute the local mass flow
394  ! for a family or sliding mesh interface.
395  ! Dimension is (1:il,2:jl,2:kl).
396  ! indFamilyJ(:,:,:) - Idem for the j-faces.
397  ! Dimension is (2:il,1:jl,2:kl).
398  ! indFamilyK(:,:,:) - Idem for the k-faces.
399  ! Dimension is (2:il,2:jl,1:kl)
400  ! factFamilyI(:,:,:) - Corresponding factor to make sure
401  ! that the massflow is defined positive
402  ! when it enters the block and to define
403  ! the mass flow of the entire wheel
404  ! instead of a sector. Hence the possible
405  ! values or -nSlices and nSlices, where
406  ! nSlices or the number of sections to
407  ! obtain the full wheel.
408  ! factFamilyJ(:,:,:) - Idem for the j-faces.
409  ! factFamilyK(:,:,:) - Idem for the k-faces.
410  !
411  ! rotMatrixI(:,:,:,:,:) - Rotation matrix of the i-faces to
412  ! transform the velocity components
413  ! from Cartesian to local cylindrical.
414  ! This is needed only for problems with
415  ! rotational periodicity in combination
416  ! with an upwind scheme.
417  ! Dimension is (1:il,2:jl,2:kl,3,3).
418  ! rotMatrixJ(:,:,:,:,:) - Idem for the j-faces.
419  ! Dimension is (2:il,1:jl,2:kl,3,3).
420  ! rotMatrixK(:,:,:,:,:) - Idem for the k-faces.
421  ! Dimension is (2:il,2:jl,1:kl,3,3).
422  !
423  ! blockIsMoving - Whether or not the block is moving.
424  ! addGridVelocities - Whether or not the face velocities
425  ! are allocated and set.
426  ! sFaceI(0:ie,je,ke) - Dot product of the face velocity and
427  ! the normal in i-direction.
428  ! sFaceJ(ie,0:je,ke) - Idem in j-direction.
429  ! sFaceK(ie,je,0:ke) - Idem in k-direction.
430 
431  real(kind=realtype), dimension(:, :, :, :), pointer :: x, xtmp
432  real(kind=realtype), dimension(:, :, :, :, :), pointer :: xold
433  real(kind=realtype), dimension(:, :, :, :), pointer :: si, sj, sk
434  real(kind=realtype), dimension(:, :, :), pointer :: vol
435  real(kind=realtype), dimension(:, :, :, :), pointer :: volold
436  real(kind=realtype), dimension(:, :, :), pointer :: volref
437  real(kind=realtype), dimension(:, :, :), pointer :: skew
438  real(kind=realtype), dimension(:, :, :, :), pointer :: uv
439  integer(kind=intType), dimension(:, :, :, :), pointer :: surfnodeindices
440 
441  integer(kind=porType), dimension(:, :, :), pointer :: pori, porj, pork
442  integer(kind=intType), dimension(:, :, :), pointer :: indfamilyi, indfamilyj, indfamilyk
443  integer(kind=intType), dimension(:, :, :), pointer :: factfamilyi, factfamilyj, factfamilyk
444  real(kind=realtype), dimension(:, :, :, :, :), pointer :: rotmatrixi, rotmatrixj, rotmatrixk
445 
446  logical :: blockismoving, addgridvelocities
447 
448  real(kind=realtype), dimension(:, :, :), pointer :: sfacei, sfacej, sfacek
449 
450  ! Added by HDN
451  ! xALE(0:ie,0:je,0:ke,3) - Temporary storage of x so that
452  ! intermediate meshes can be stored in
453  ! x directly
454  ! sVeloIALE(0:ie,1:je,1:ke,3)
455  ! sVeloJALE(1:ie,0:je,1:ke,3)
456  ! sVeloKALE(1:ie,1:je,0:ke,3) - Storage of surface velocities at one
457  ! time step
458  ! sIALE(0:nALEsteps,0:ie,1:je,1:ke,3)
459  ! sJALE(0:nALEsteps,1:ie,0:je,1:ke,3)
460  ! sKALE(0:nALEsteps,1:ie,1:je,0:ke,3) - Storage of sI, sJ, sK for intermediate
461  ! meshes
462  ! sFaceIALE(0:nALEsteps,0:ie,je,ke)
463  ! sFaceJALE(0:nALEsteps,ie,0:je,ke)
464  ! sFaceKALE(0:nALEsteps,ie,je,0:ke) - Storage of sFaceI, sFaceJ, sFaceK for
465  ! intermediate meshes
466  real(kind=realtype), dimension(:, :, :, :), pointer :: xale
467  real(kind=realtype), dimension(:, :, :, :), pointer :: sveloiale, svelojale, svelokale
468  real(kind=realtype), dimension(:, :, :, :, :), pointer :: siale, sjale, skale
469  real(kind=realtype), dimension(:, :, :, :), pointer :: sfaceiale, sfacejale, sfacekale
470 
471  ! Tempory storage for overset variables
472  real(kind=realtype), dimension(:, :, :, :), pointer :: xseed
473  integer(kind=intType), dimension(:, :, :), pointer :: wallind
474 
475  !
476  ! Flow variables.
477  !
478  ! w(0:ib,0:jb,0:kb,1:nw) - The set of independent variables
479  ! w(i,j,k,1:nwf) flow field
480  ! variables, which are rho, u,
481  ! v, w and rhoE. In other words
482  ! the velocities are stored and
483  ! not the momentum!!!!
484  ! w(i,j,k,nt1:nt2) turbulent
485  ! variables; also the primitive
486  ! variables are stored.
487  ! wOld(nOld,2:il,2:jl,2:kl,nw) - Solution on older time levels,
488  ! needed for the time integration
489  ! for unsteady problems. In
490  ! constrast to w, the conservative
491  ! variables are stored in wOld for
492  ! the flow variables; the turbulent
493  ! variables are always the
494  ! primitive ones.
495  ! Only allocated on the finest
496  ! mesh.
497  ! p(0:ib,0:jb,0:kb) - Static pressure.
498  ! gamma(0:ib,0:jb,0:kb) - Specific heat ratio; only
499  ! allocated on the finest grid.
500  ! rlv(0:ib,0:jb,0:kb) - Laminar viscosity; only
501  ! allocated on the finest mesh
502  ! and only for viscous problems.
503  ! rev(0:ib,0:jb,0:kb) - Eddy viscosity; only
504  ! allocated rans problems with
505  ! eddy viscosity models.
506  ! s(1:ie,1:je,1:ke,3) - Mesh velocities of the cell
507  ! centers; only for moving mesh
508  ! problems.
509 
510  real(kind=realtype), dimension(:, :, :, :), pointer :: w, wtmp
511  real(kind=realtype), dimension(:, :, :, :, :), pointer :: dw_deriv
512  real(kind=realtype), dimension(:, :, :, :, :), pointer :: wold
513  real(kind=realtype), dimension(:, :, :), pointer :: p, gamma, aa
514  real(kind=realtype), dimension(:, :, :), pointer :: rlv, rev
515  real(kind=realtype), dimension(:, :, :, :), pointer :: s
516  real(kind=realtype), dimension(:, :, :), pointer :: shocksensor
517 
518  ! Nodal Fluxes: ux,uy,uz,vx,vy,vz,wx,wy,wz,qx,qy,qz(il, jl, kl)
519  real(kind=realtype), dimension(:, :, :), pointer :: ux, uy, uz
520  real(kind=realtype), dimension(:, :, :), pointer :: vx, vy, vz
521  real(kind=realtype), dimension(:, :, :), pointer :: wx, wy, wz
522  real(kind=realtype), dimension(:, :, :), pointer :: qx, qy, qz
523 
524  !
525  ! Residual and multigrid variables.
526  !
527  ! dw(0:ib,0:jb,0:kb,1:nw) - Values of convective and combined
528  ! flow residuals. Only allocated on
529  ! the finest mesh.
530  ! fw(0:ib,0:jb,0:kb,1:nwf) - values of artificial dissipation
531  ! and viscous residuals.
532  ! Only allocated on the finest mesh.
533 
534  ! dwOldRK(:,2:il,2:jl,2:kl,nw) - Old residuals for the time
535  ! accurate Runge-Kutta schemes.
536  ! The first dimension is
537  ! nRKStagesUnsteady - 1.Only
538  ! allocated on the finest level
539  ! and only in unsteady mode for
540  ! Runge-Kutta schemes.
541 
542  ! w1(1:ie,1:je,1:ke,1:nMGVar) - Values of the mg variables
543  ! upon first entry to a coarser
544  ! mesh; only allocated on the
545  ! coarser grids. The variables
546  ! used to compute the multigrid
547  ! corrections are rho, u, v, w
548  ! and p; the rhoE value is used
549  ! for unsteady problems only.
550  ! p1(1:ie,1:je,1:ke) - Value of the pressure upon
551  ! first entry to a coarser grid;
552  ! only allocated on the coarser
553  ! grids.
554  ! wr(2:il,2:jl,2:kl,1:nMGVar) - Multigrid forcing terms; only
555  ! allocated on the coarser grids.
556  ! The forcing term of course
557  ! contains conservative residuals,
558  ! at least for the flow variables.
559  ! shockSensor(0:ib,0:jb,0:kb) Precomputed sensor value for shock
560  ! that is *NOT* differentated.
561  ! scratch(0:ib,0:jb,0:kb,5) Scratch space for the turbulence
562  ! models. NOMINALLY this could use
563  ! dw and the code was nominally setup
564  ! for this originally. However, this
565  ! complicates reverse mode sensitivities
566  ! So we use this instead.
567 
568  real(kind=realtype), dimension(:, :, :), pointer :: p1
569  real(kind=realtype), dimension(:, :, :, :), pointer :: dw, fw
570  real(kind=realtype), dimension(:, :, :, :), pointer :: dwtmp, dwtmp2
571  real(kind=realtype), dimension(:, :, :, :, :), pointer :: dwoldrk
572  real(kind=realtype), dimension(:, :, :, :), pointer :: w1, wr
573  real(kind=realtype), dimension(:, :, :, :), pointer :: scratch
574 
575  ! Added by HDN
576  ! Used for ALE. Only allocated on the finest mesh.
577  ! Extra dim is used to store initial residuals
578  ! dwALE(0:nALEsteps,0:ib,0:jb,0:kb,1:nw) - Values of ONLY the convective flux
579  ! of intermediate meshes.
580  ! fwALE(0:nALEsteps,0:ib,0:jb,0:kb,1:nwf) - values of ONLY the artificial
581  ! dissipation of intermediate meshes.
582  real(kind=realtype), dimension(:, :, :, :, :), pointer :: dwale, fwale
583 
584  ! mgIFine(2:il,2) - The two fine grid i-cells used for the
585  ! restriction of the solution and residual to
586  ! the coarse grid. Only on the coarser grids.
587  ! mgJFine(2:jl,2) - Idem for j-cells.
588  ! mgKFine(2:kl,2) - Idem for k-cells.
589 
590  ! mgIWeight(2:il) - Weight for the residual restriction in
591  ! in i-direction. Value is either 0.5 or 1.0,
592  ! depending whether mgIFine(,1) is equal to
593  ! or differs from mgIFine(,2).
594  ! mgJWeight(2:jl) - Idem for weights in j-direction.
595  ! mgKWeight(2:kl) - Idem for weights in k-direction.
596 
597  ! mgICoarse(2:il,2) - The two coarse grid i-cells used for the
598  ! interpolation of the correction to the
599  ! fine grid. Not on the coarsest grid.
600  ! mgJCoarse(2:jl,2) - Idem for j-cells.
601  ! mgKCoarse(2:kl,2) - Idem for k-cells.
602 
603  integer(kind=intType), dimension(:, :), pointer :: mgifine
604  integer(kind=intType), dimension(:, :), pointer :: mgjfine
605  integer(kind=intType), dimension(:, :), pointer :: mgkfine
606 
607  real(kind=realtype), dimension(:), pointer :: mgiweight
608  real(kind=realtype), dimension(:), pointer :: mgjweight
609  real(kind=realtype), dimension(:), pointer :: mgkweight
610 
611  integer(kind=intType), dimension(:, :), pointer :: mgicoarse
612  integer(kind=intType), dimension(:, :), pointer :: mgjcoarse
613  integer(kind=intType), dimension(:, :), pointer :: mgkcoarse
614 
615  ! iCoarsened - How this block was coarsened in i-direction.
616  ! jCoarsened - How this block was coarsened in j-direction.
617  ! kCoarsened - How this block was coarsened in k-direction.
618 
619  integer(kind=porType) :: icoarsened, jcoarsened, kcoarsened
620 
621  ! iCo: Indicates whether or not i grid lines are present on the
622  ! coarse grid; not allocated for the coarsest grid.
623  ! jCo: Idem in j-direction.
624  ! kCo: Idem in k-direction.
625 
626  logical, dimension(:), pointer :: ico, jco, kco
627  !
628  ! Time-stepping and spectral radii variables.
629  ! only allocated on the finest grid.
630  !
631  ! wn(2:il,2:jl,2:kl,1:nMGVar) - Values of the update variables
632  ! at the beginning of the RungeKutta
633  ! iteration. Only allocated for
634  ! RungeKutta smoother.
635  ! pn(2:il,2:jl,2:kl) - The pressure for the RungeKutta
636  ! smoother.
637  ! dtl(1:ie,1:je,1:ke) - Time step
638  ! radI(1:ie,1:je,1:ke) - Spectral radius in i-direction.
639  ! radJ(1:ie,1:je,1:ke) - Spectral radius in j-direction.
640  ! radK(1:ie,1:je,1:ke) - Spectral radius in k-direction.
641 
642  real(kind=realtype), dimension(:, :, :, :), pointer :: wn
643  real(kind=realtype), dimension(:, :, :), pointer :: pn
644  real(kind=realtype), dimension(:, :, :), pointer :: dtl
645  real(kind=realtype), dimension(:, :, :), pointer :: radi, radj, radk
646 
647  !
648  ! Variables for Iso/Surface Slice generation
649  ! fc(1:ie,1:je,1:ke) - cell center values of the function to be iso-valued
650  ! fn(1:il,1:jl,1:kl) - node values of the function to be iso-valued
651  ! Note these are are only allocated temporaily during solution writing.
652 
653  real(kind=realtype), dimension(:, :, :), pointer :: fc
654  real(kind=realtype), dimension(:, :, :), pointer :: fn
655 
656  !
657  ! Turbulence model variables.
658  !
659  ! d2Wall(2:il,2:jl,2:kl) - Distance from the center of the cell
660  ! to the nearest viscous wall.
661  ! intermittency( ) - Function defining the transition location
662 
663  real(kind=realtype), dimension(:, :, :), pointer :: d2wall, filterdes
664  real(kind=realtype), dimension(:, :, :), pointer :: intermittency
665 
666  ! bmti1(je,ke,nt1:nt2,nt1:nt2): Matrix used for the implicit
667  ! boundary condition treatment of
668  ! the turbulence equations at the
669  ! iMin boundary. Only allocated on
670  ! the finest level and for the 1st
671  ! spectral solution.
672  ! bmti2(je,ke,nt1:nt2,nt1:nt2): Idem for the iMax boundary.
673  ! bmtj1(ie,ke,nt1:nt2,nt1:nt2): Idem for the jMin boundary.
674  ! bmtj2(ie,ke,nt1:nt2,nt1:nt2): Idem for the jMax boundary.
675  ! bmtk1(ie,je,nt1:nt2,nt1:nt2): Idem for the kMin boundary.
676  ! bmtk2(ie,je,nt1:nt2,nt1:nt2): Idem for the kMax boundary.
677 
678  real(kind=realtype), dimension(:, :, :, :), pointer :: bmti1, bmti2
679  real(kind=realtype), dimension(:, :, :, :), pointer :: bmtj1, bmtj2
680  real(kind=realtype), dimension(:, :, :, :), pointer :: bmtk1, bmtk2
681 
682  ! bvti1(je,ke,nt1:nt2): RHS vector used for the implicit
683  ! boundary condition treatment of the
684  ! turbulence equations at the iMin
685  ! boundary. Only allocated on the finest
686  ! level and for the 1st spectral solution.
687  ! bvti2(je,ke,nt1:nt2): Idem for the iMax boundary.
688  ! bvtj1(ie,ke,nt1:nt2): Idem for the jMin boundary.
689  ! bvtj2(ie,ke,nt1:nt2): Idem for the jMax boundary.
690  ! bvti2(je,ke,nt1:nt2): Idem for the iMax boundary.
691  ! bvtk1(ie,ke,nt1:nt2): Idem for the kMin boundary.
692  ! bvtk2(ie,ke,nt1:nt2): idem for the kMax boundary.
693 
694  real(kind=realtype), dimension(:, :, :), pointer :: bvti1, bvti2
695  real(kind=realtype), dimension(:, :, :), pointer :: bvtj1, bvtj2
696  real(kind=realtype), dimension(:, :, :), pointer :: bvtk1, bvtk2
697  !
698  ! Relation to the original cgns grid.
699  !
700  ! sectionID - The section of the grid this block belongs to.
701  ! cgnsBlockID - Block/zone number of the cgns grid to which
702  ! this block is related.
703  ! iBegOr, iEndOr - Range of points of this block in the
704  ! jBegOr, jEndOr corresponding cgns block, i.e. for this block
705  ! kBegOr, kEndOr iBegOr <= i <= iEndOr, jBegOr <= j <= jEndOr,
706  ! kBegOr <= k <= kEndOr.
707  ! It is of course possible that the entire
708  ! block is stored.
709  integer(kind=intType) :: sectionid = 1
710  integer(kind=intType) :: cgnsblockid
711  integer(kind=intType) :: ibegor, iendor, jbegor, jendor
712  integer(kind=intType) :: kbegor, kendor
713  type(surfacenodeweightarray), dimension(6) :: nodalweights
714  !
715  ! Adjoint solver variables.
716  !
717  ! globalNode(ib:ie,jb:je,kb:ke): Global node numbering.
718  ! globalCell(0:ib,0:jb,0:kb): Global cell numbering.
719  ! color(0:ib,0:jb,0:kb) : Temporary coloring array used for
720  ! forward mode AD/FD calculations
721  integer(kind=intType), dimension(:, :, :), pointer :: globalcgnsnode
722  integer(kind=intType), dimension(:, :, :), pointer :: globalnode
723  integer(kind=intType), dimension(:, :, :), pointer :: globalcell
724  integer(kind=intType), dimension(:, :, :), pointer :: color
725 
726  ! Data storing the first order PC in tri-diagonal ordering. 7
727  ! real(kind=realType), dimension(:, :, :, :, :), pointer :: Diag
728  ! real(kind=realType), dimension(:, :, :, :, :), pointer :: i_L, i_U
729  ! real(kind=realType), dimension(:, :, :, :, :), pointer :: j_L, j_U
730  ! real(kind=realType), dimension(:, :, :, :, :), pointer :: k_L, k_U
731 
732  real(kind=realtype), dimension(:, :, :, :, :), pointer :: pcmat
733 
734  ! Generic vectors for doing products/preconditioning. Like w, but
735  ! only 1 level of halos, and it is in block ordering (nw first)
736  ! instead of field ordering like w is.
737  real(kind=realtype), dimension(:, :, :, :), pointer :: pcvec1, pcvec2
738 
739  ! Data for the factorized trigonal solves
740  real(kind=realtype), dimension(:, :, :, :), pointer :: i_d_fact, j_d_fact, k_d_fact
741  real(kind=realtype), dimension(:, :, :, :), pointer :: i_l_fact, j_l_fact, k_l_fact
742  real(kind=realtype), dimension(:, :, :, :), pointer :: i_u_fact, j_u_fact, k_u_fact
743  real(kind=realtype), dimension(:, :, :, :), pointer :: i_u2_fact, j_u2_fact, k_u2_fact
744 
745  integer(kind=intType), dimension(:, :, :, :), pointer :: i_ipiv, j_ipiv, k_ipiv
746 
747  ! A list of pointers for generic communication of either real or
748  ! integer data.
749  type(rptr), dimension(24) :: realcommvars
750  type(iptr), dimension(3) :: intcommvars
751 
752  end type blocktype
753 
754  !
755  ! Array of all blocks at all multigrid levels and spectral sols.
756  !
757  ! nDom: total number of computational blocks.
758  ! flowDoms(:,:,:): array of blocks. Dimensions are
759  ! (nDom,nLevels,nTimeIntervalsSpectral)
760 
761  integer(kind=intType) :: ndom
762 
763  ! A global paramter for how to sort fringes
764  integer(kind=intType) :: fringesorttype = sortbydonor
765 
766 #ifdef USE_TAPENADE
767  ! This is never actually compiled...just make tapenade think it
768  ! isn't allocatable
769  type(blocktype), dimension(nn:nn, 1, ntimeIntervalsSpectral) :: flowdoms
770 #else
771  type(blocktype), allocatable, target, dimension(:, :, :) :: flowdoms
772  type(blocktype), allocatable, target, dimension(:, :, :) :: flowdomsd
773  type(blocktype), allocatable, target, dimension(:, :, :) :: flowdomsb
774 #endif
775 
776  !
777  ! Additional info needed in the flow solver.
778  !
779  ! nCellGlobal(nLev) - Global number of cells on every mg level.
780 
781  integer(kind=intType), allocatable, dimension(:) :: ncellglobal
782 
783 contains
784 
785  logical function lessequalfringetype(g1, g2)
786 
787  ! lessEqual returns .true. if g1 <= g2 and .false. otherwise.
788  ! The comparison is firstly based on the processor ID of the
789  ! donor, then the block, then then the I, J, K
790  !
791  implicit none
792  !
793  ! Function arguments.
794  !
795  type(fringetype), intent(in) :: g1, g2
796  !
797  ! Compare the donor processors first. If not equal,
798  ! set lessEqual appropriately and return.
799  if (fringesorttype == sortbydonor) then
800  if (g1%donorProc < g2%donorProc) then
801  lessequalfringetype = .true.
802  return
803  else if (g1%donorProc > g2%donorProc) then
804  lessequalfringetype = .false.
805  return
806  end if
807 
808  ! Donor processors are identical. Now we check the block
809 
810  if (g1%donorBlock < g2%donorBlock) then
811  lessequalfringetype = .true.
812  return
813  else if (g1%donorBlock > g2%donorBlock) then
814  lessequalfringetype = .false.
815  return
816  end if
817 
818  ! Compare the indices of the halo. First k, then j and
819  ! finally i.
820 
821  if (g1%dIndex < g2%dIndex) then
822  lessequalfringetype = .true.
823  return
824  else if (g1%dindex > g2%dIndex) then
825  lessequalfringetype = .false.
826  return
827  end if
828 
829  else if (fringesorttype == sortbyreceiver) then
830 
831  ! Compare my indices
832 
833  if (g1%myIndex < g2%myIndex) then
834  lessequalfringetype = .true.
835  return
836  else if (g1%myIndex > g2%myIndex) then
837  lessequalfringetype = .false.
838  return
839  end if
840 
841  ! Now compare the donor information:
842 
843  if (g1%donorProc < g2%donorProc) then
844  lessequalfringetype = .true.
845  return
846  else if (g1%donorProc > g2%donorProc) then
847  lessequalfringetype = .false.
848  return
849  end if
850 
851  ! Donor processors are identical. Now we check the block
852 
853  if (g1%donorBlock < g2%donorBlock) then
854  lessequalfringetype = .true.
855  return
856  else if (g1%donorBlock > g2%donorBlock) then
857  lessequalfringetype = .false.
858  return
859  end if
860 
861  ! Compare the indices of the halo. First k, then j and
862  ! finally i.
863 
864  if (g1%dIndex < g2%dIndex) then
865  lessequalfringetype = .true.
866  return
867  else if (g1%dIndex > g2%dIndex) then
868  lessequalfringetype = .false.
869  return
870  end if
871  end if
872 
873  ! Both entities are identical. So set lessEqual to .true.
874 
875  lessequalfringetype = .true.
876 
877  end function lessequalfringetype
878 
879  logical function lessfringetype(g1, g2)
880 
881  ! less returns .true. if g1 <= g2 and .false. otherwise.
882  ! The comparison is firstly based on the processor ID of the
883  ! donor, then the block, then then the I, J, K
884  !
885  implicit none
886  !
887  ! Function arguments.
888  !
889  type(fringetype), intent(in) :: g1, g2
890  !
891  ! Compare the donor processors first. If not equal,
892  ! set less appropriately and return.
893  if (fringesorttype == sortbydonor) then
894  if (g1%donorProc < g2%donorProc) then
895  lessfringetype = .true.
896  return
897  else if (g1%donorProc > g2%donorProc) then
898  lessfringetype = .false.
899  return
900  end if
901 
902  ! Donor processors are identical. Now we check the block
903 
904  if (g1%donorBlock < g2%donorBlock) then
905  lessfringetype = .true.
906  return
907  else if (g1%donorBlock > g2%donorBlock) then
908  lessfringetype = .false.
909  return
910  end if
911 
912  ! Compare the indices of the halo. First k, then j and
913  ! finally i.
914 
915  if (g1%dIndex < g2%dIndex) then
916  lessfringetype = .true.
917  return
918  else if (g1%dIndex > g2%dIndex) then
919  lessfringetype = .false.
920  return
921  end if
922 
923  else if (fringesorttype == sortbyreceiver) then
924 
925  ! Compare my indices
926 
927  if (g1%myIndex < g2%myIndex) then
928  lessfringetype = .true.
929  return
930  else if (g1%myIndex > g2%myIndex) then
931  lessfringetype = .false.
932  return
933  end if
934 
935  ! Now compare the donor information:
936 
937  if (g1%donorProc < g2%donorProc) then
938  lessfringetype = .true.
939  return
940  else if (g1%donorProc > g2%donorProc) then
941  lessfringetype = .false.
942  return
943  end if
944 
945  ! Donor processors are identical. Now we check the block
946 
947  if (g1%donorBlock < g2%donorBlock) then
948  lessfringetype = .true.
949  return
950  else if (g1%donorBlock > g2%donorBlock) then
951  lessfringetype = .false.
952  return
953  end if
954 
955  ! Compare the indices of the halo. First k, then j and
956  ! finally i.
957 
958  if (g1%dIndex < g2%dIndex) then
959  lessfringetype = .true.
960  return
961  else if (g1%dIndex > g2%dIndex) then
962  lessfringetype = .false.
963  return
964  end if
965  end if
966 
967  ! Both entities are identical. So set less to .False.
968 
969  lessfringetype = .false.
970 
971  end function lessfringetype
972 
973  logical function lessequalinterppttype(g1, g2)
974 
975  ! lessEqual returns .true. if g1 <= g2 and .false. otherwise.
976  ! The comparison is firstly based on the processor ID of the
977  ! donor, then the block, then then the I, J, K
978  !
979  implicit none
980  !
981  ! Function arguments.
982  !
983  type(interppttype), intent(in) :: g1, g2
984  !
985 
986  if (g1%donorProc < g2%donorProc) then
987  lessequalinterppttype = .true.
988  return
989  else if (g1%donorProc > g2%donorProc) then
990  lessequalinterppttype = .false.
991  return
992  end if
993 
994  ! Donor processors are identical. Now we check the block
995 
996  if (g1%donorBlock < g2%donorBlock) then
997  lessequalinterppttype = .true.
998  return
999  else if (g1%donorBlock > g2%donorBlock) then
1000  lessequalinterppttype = .false.
1001  return
1002  end if
1003 
1004  ! Compare the indices of the halo. First k, then j and
1005  ! finally i.
1006 
1007  if (g1%dK < g2%dK) then
1008  lessequalinterppttype = .true.
1009  return
1010  else if (g1%dK > g2%dK) then
1011  lessequalinterppttype = .false.
1012  return
1013  end if
1014 
1015  if (g1%dJ < g2%dJ) then
1016  lessequalinterppttype = .true.
1017  return
1018  else if (g1%dJ > g2%dJ) then
1019  lessequalinterppttype = .false.
1020  return
1021  end if
1022 
1023  if (g1%dI < g2%dI) then
1024  lessequalinterppttype = .true.
1025  return
1026  else if (g1%dI > g2%dI) then
1027  lessequalinterppttype = .false.
1028  return
1029  end if
1030 
1031  ! Both entities are identical. So set lessEqual to .true.
1032 
1033  lessequalinterppttype = .true.
1034 
1035  end function lessequalinterppttype
1036 
1037  logical function lessinterppttype(g1, g2)
1038 
1039  implicit none
1040  !
1041  ! Function arguments.
1042  !
1043  type(interppttype), intent(in) :: g1, g2
1044  !
1045  if (g1%donorProc < g2%donorProc) then
1046  lessinterppttype = .true.
1047  return
1048  else if (g1%donorProc > g2%donorProc) then
1049  lessinterppttype = .false.
1050  return
1051  end if
1052 
1053  ! Donor processors are identical. Now we check the block
1054 
1055  if (g1%donorBlock < g2%donorBlock) then
1056  lessinterppttype = .true.
1057  return
1058  else if (g1%donorBlock > g2%donorBlock) then
1059  lessinterppttype = .false.
1060  return
1061  end if
1062 
1063  ! Compare the indices of the halo. First k, then j and
1064  ! finally i.
1065 
1066  if (g1%dK < g2%dK) then
1067  lessinterppttype = .true.
1068  return
1069  else if (g1%dK > g2%dK) then
1070  lessinterppttype = .false.
1071  return
1072  end if
1073 
1074  if (g1%dJ < g2%dJ) then
1075  lessinterppttype = .true.
1076  return
1077  else if (g1%dJ > g2%dJ) then
1078  lessinterppttype = .false.
1079  return
1080  end if
1081 
1082  if (g1%dI < g2%dI) then
1083  lessinterppttype = .true.
1084  return
1085  else if (g1%dI > g2%dI) then
1086  lessinterppttype = .false.
1087  return
1088  end if
1089 
1090  ! Both entities are identical. So set less to .False.
1091 
1092  lessinterppttype = .false.
1093 
1094  end function lessinterppttype
1095 
1096 end module block
Definition: BCData.F90:1
Definition: block.F90:1
integer(kind=inttype) ndom
Definition: block.F90:761
logical function lessequalinterppttype(g1, g2)
Definition: block.F90:974
type(blocktype), dimension(:, :, :), allocatable, target flowdomsb
Definition: block.F90:773
type(blocktype), dimension(:, :, :), allocatable, target flowdomsd
Definition: block.F90:772
logical function lessequalfringetype(g1, g2)
Definition: block.F90:786
integer(kind=inttype) fringesorttype
Definition: block.F90:764
logical function lessinterppttype(g1, g2)
Definition: block.F90:1038
logical function lessfringetype(g1, g2)
Definition: block.F90:880
type(blocktype), dimension(:, :, :), allocatable, target flowdoms
Definition: block.F90:771
integer(kind=inttype), dimension(:), allocatable ncellglobal
Definition: block.F90:781
integer(kind=inttype), parameter sortbyreceiver
Definition: constants.F90:325
integer, parameter maxcgnsnamelen
Definition: constants.F90:17
integer(kind=inttype), parameter sortbydonor
Definition: constants.F90:324