FLASH-X
Doxygen Generated Documentation From Interface Source Code
Grid_interface.F90
Go to the documentation of this file.
1!! NOTICE
2!! Copyright 2022 UChicago Argonne, LLC and contributors
3!!
4!! Licensed under the Apache License, Version 2.0 (the "License");
5!! you may not use this file except in compliance with the License.
6!!
7!! Unless required by applicable law or agreed to in writing, software
8!! distributed under the License is distributed on an "AS IS" BASIS,
9!! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
10!! See the License for the specific language governing permissions and
11!! limitations under the License.
12!!
13!! This is the header file for the Grid unit that defines its
14!! public interfaces.
15
16!!REORDER(4): fluxBuf[XYZ], fluxCorr[XYZ], fluxOld[XYZ]
17
19
23
24 implicit none
25
26#include "Simulation.h"
27#ifdef Grid_releaseBlkPtr
28! disabling drift macro expansion because it doesn't apply here, only to call sites. see: drift
29#undef Grid_releaseBlkPtr
30#endif
31
32#include "constants.h"
33#include "Particles.h"
34!#include "GridParticles.h"
35
36 integer,parameter :: GRID_PDE_BND_ISOLATED = 0
37 integer,parameter :: GRID_PDE_BND_PERIODIC = 1
38 integer,parameter :: GRID_PDE_BND_DIRICHLET = 2
39 integer,parameter :: GRID_PDE_BND_NEUMANN = 3
40 integer,parameter :: GRID_PDE_BND_GIVENVAL = 4
41 integer,parameter :: GRID_PDE_BND_GIVENGRAD = 5
42
43 integer,parameter :: GRID_COPYDIR_TO_VECT = 1
44 integer,parameter :: GRID_COPYDIR_FROM_VECT = 2
45
46#include "FortranLangFeatures.fh"
47
49 subroutine Grid_ascGetBlkPtr(blockID,dataPtr, gridDataStruct)
50 implicit none
51 integer, intent(in) :: blockID
52 real, dimension(:,:,:,:), pointer :: dataPtr
53 integer, optional,intent(in) :: gridDataStruct
54 end subroutine Grid_ascGetBlkPtr
55 subroutine Grid_ascGetBlk5Ptr(blockID,data5Ptr, gridDataStruct)
56 implicit none
57 integer, intent(in) :: blockID
58 real, dimension(:,:,:,:,:), pointer :: data5Ptr
59 integer, optional,intent(in) :: gridDataStruct
60 end subroutine Grid_ascGetBlk5Ptr
61 end interface Grid_ascGetBlkPtr
62
63 interface
64 subroutine Grid_ascReleaseBlkPtr(blockId, dataPtr, gridDataStruct)
65 implicit none
66 integer,intent(in) :: blockId
67 real, pointer :: dataPtr(:,:,:,:)
68 integer,optional, intent(in) :: gridDataStruct
69 end subroutine Grid_ascReleaseBlkPtr
70 subroutine Grid_ascReleaseBlk5Ptr(blockId, data5Ptr, gridDataStruct)
71 implicit none
72 integer,intent(in) :: blockId
73 real, POINTER_INTENT_OUT :: data5Ptr(:,:,:,:,:)
74 integer,optional, intent(in) :: gridDataStruct
75 end subroutine Grid_ascReleaseBlk5Ptr
76 end interface
77
78 interface
79 subroutine Grid_genGetBlkPtr(blockID,dataPtr, varDesc,dataPtr2,dataPtr3)
80 implicit none
81 integer, intent(in) :: blockID
82 real, dimension(:,:,:,:), POINTER_INTENT_OUT :: dataPtr
83 integer,intent(in),OPTIONAL :: varDesc(:)
84 real, dimension(:,:,:,:), POINTER_INTENT_OUT,OPTIONAL :: dataPtr2, dataPtr3
85 end subroutine Grid_genGetBlkPtr
86 end interface
87
88 interface
89 subroutine Grid_genReleaseBlkPtr(blockID, dataPtr, varDesc,dataPtr2,dataPtr3)
90 implicit none
91 integer,intent(in) :: blockId
92 real, pointer :: dataPtr(:,:,:,:)
93 integer,intent(in),OPTIONAL :: varDesc(:)
94 real, dimension(:,:,:,:), pointer,OPTIONAL :: dataPtr2, dataPtr3
95 end subroutine Grid_genReleaseBlkPtr
96 end interface
97
98 interface
99 subroutine Grid_applyBCEdge(bcType,bcDir,guard,var,dataRow,face,&
100 gridDataStruct, secondCoord, thirdCoord)
101 integer,intent(IN):: bcType,bcDir,guard,var,face,gridDataStruct
102 real,dimension(:),intent(INOUT)::dataRow
103 real,intent(IN),OPTIONAL :: secondCoord,thirdCoord
104 end subroutine Grid_applyBCEdge
105 end interface
106
107 interface
108 subroutine Grid_applyBCEdgeAllUnkVars(bcType,bcDir,guard,dataRow,face,&
109 cellCenterSweepCoord, secondCoord,thirdCoord)
110 integer,intent(IN):: bcType,bcDir,guard,face
111 real,dimension(2*guard,NUNK_VARS),intent(INOUT)::dataRow
112 real,intent(IN):: cellCenterSweepCoord(*), secondCoord,thirdCoord
113 end subroutine Grid_applyBCEdgeAllUnkVars
114 end interface
115
118 end subroutine Grid_computeUserVars
119 end interface
120
121 interface
122 subroutine Grid_computeVarMean(iUnk, mean)
123 implicit none
124 integer, intent(in) :: iUnk
125 real, intent(out) :: mean
126 end subroutine Grid_computeVarMean
127 end interface
128
129 interface
130 subroutine Grid_computeVarNorm (level, normType, ivar, norm, leafOnly)
131 implicit none
132 integer, intent(IN) :: normType, level, ivar, leafOnly
133 real, intent(OUT) :: norm
134 end subroutine Grid_computeVarNorm
135 end interface
136
137 interface
138 subroutine Grid_computeVarDiff(level, gr_iRefSoln, gr_iSoln, ires)
139 implicit none
140 integer, intent(in) :: level, gr_iRefSoln, gr_iSoln, ires
141 end subroutine Grid_computeVarDiff
142 end interface
143
144 interface
145 subroutine Grid_communicateFluxes(axis, coarse_level)
146 implicit none
147 integer, intent(IN) :: axis
148 integer, intent(IN) :: coarse_level
149 end subroutine Grid_communicateFluxes
150 end interface
151
153 subroutine Grid_conserveFluxes(axis, coarse_level)
154 integer, intent(IN) :: axis
155 integer, intent(IN) :: coarse_level
156 end subroutine Grid_conserveFluxes
157 end interface
158
159 interface Grid_dump
160 subroutine Grid_dump(var,num,solnData,blockDesc,gcell)
161 use Grid_tile, ONLY : Grid_tile_t
162 implicit none
163 integer, intent(IN) :: num
164 integer, dimension(num), intent(IN) :: var
165 real,dimension(:,:,:,:),pointer :: solnData
166 type(Grid_tile_t),target, intent(in) :: blockDesc
167 logical, intent(IN) :: gcell
168 end subroutine Grid_dump
169 end interface
170
171
172 interface
173 subroutine Grid_fillGuardCells( gridDataStruct,idir,&
174 minLayers,eosMode,doEos,maskSize,mask,makeMaskConsistent,&
175 doLogMask,selectBlockType,unitReadsMeshDataOnly)
176 integer, intent(in) :: gridDataStruct
177 integer, intent(in) :: idir
178 integer,optional,intent(in) :: minLayers
179 integer,optional,intent(in) :: eosMode
180 logical,optional,intent(IN) :: doEos
181 integer,optional,intent(in) :: maskSize
182 logical,optional,dimension(:), intent(IN) :: mask
183 logical,optional, intent(IN) :: makeMaskConsistent
184 logical,optional,intent(IN) :: doLogMask
185 integer,optional,intent(in) :: selectBlockType
186 logical, optional, intent(IN) :: unitReadsMeshDataOnly
187 end subroutine Grid_fillGuardCells
188 end interface
189
191 subroutine Grid_finalize()
192 end subroutine Grid_finalize
193 end interface
194
196 subroutine Grid_getBlkCenterCoords(blockDesc,blockCenter)
197 use Grid_tile, ONLY : Grid_tile_t
198 implicit none
199 type(Grid_tile_t), intent(in) :: blockDesc
200 real, dimension(MDIM), intent(out) :: blockCenter
201 end subroutine Grid_getBlkCenterCoords
202 end interface
203
205 subroutine Grid_getBlkCornerID_desc(blockDesc, cornerID, stride,cornerIDHigh,inRegion)
206 use Grid_tile, ONLY : Grid_tile_t
207 implicit none
208 type(Grid_tile_t), intent(in) :: blockDesc
209 integer,dimension(MDIM), intent(OUT) :: cornerID, stride
210 integer,dimension(MDIM),optional,intent(OUT) :: cornerIDHigh
211 logical, optional, intent(IN) :: inRegion
212 end subroutine Grid_getBlkCornerID_desc
213 subroutine Grid_getBlkCornerID(blockID, cornerID, stride,cornerIDHigh,inRegion)
214 implicit none
215 integer, intent(in) :: blockID
216 integer,dimension(MDIM), intent(OUT) :: cornerID, stride
217 integer,dimension(MDIM),optional,intent(OUT) :: cornerIDHigh
218 logical, optional, intent(IN) :: inRegion
219 end subroutine Grid_getBlkCornerID
220 end interface
221
222
223 interface
224 subroutine Grid_getBlkIndexLimits_STANDALONE(blockId, blkLimits, blkLimitsGC,gridDataStruct)
225 integer,intent(IN) :: blockId
226 integer,dimension(2,MDIM), intent(OUT) :: blkLimits, blkLimitsGC
227 integer,optional,intent(IN) :: gridDataStruct
229 end interface
230
232 subroutine Grid_getCellCoords(axis, edge, level, lo, hi, coordinates)
233 integer, intent(in) :: axis
234 integer, intent(in) :: edge
235 integer, intent(in) :: level
236 integer, intent(in) :: lo(1:MDIM)
237 integer, intent(in) :: hi(1:MDIM)
238 real, intent(out) :: coordinates(:)
239 end subroutine Grid_getCellCoords
240 end interface
241
243 subroutine Grid_getCellFaceAreas(axis, level, lo, hi, areas)
244 integer, intent(in) :: axis
245 integer, intent(in) :: level
246 integer, intent(in) :: lo(1:MDIM)
247 integer, intent(in) :: hi(1:MDIM)
248 real, intent(out) :: areas(lo(IAXIS):hi(IAXIS), &
249 lo(JAXIS):hi(JAXIS), &
250 lo(KAXIS):hi(KAXIS))
251 end subroutine Grid_getCellFaceAreas
252 end interface
253
255 subroutine Grid_getCellVolumes(level, lo, hi, volumes)
256 integer, intent(in) :: level
257 integer, intent(in) :: lo(1:MDIM)
258 integer, intent(in) :: hi(1:MDIM)
259 real, intent(out) :: volumes(lo(IAXIS):hi(IAXIS), &
260 lo(JAXIS):hi(JAXIS), &
261 lo(KAXIS):hi(KAXIS))
262 end subroutine Grid_getCellVolumes
263 end interface
264
266 subroutine Grid_getDeltas(lev, del)
267 integer, intent(in) :: lev
268 real, dimension(MDIM), intent(out) :: del
269 end subroutine Grid_getDeltas
270 end interface
271
272 interface
273 subroutine Grid_getGlobalIndexLimits(globalIndexLimits)
274 integer, dimension(MDIM), intent(out) :: globalIndexLimits
275 end subroutine Grid_getGlobalIndexLimits
276 end interface
277
278 interface
279 subroutine Grid_getListOfBlocks(blockType, listOfBlocks,count,refinementLevel,region_bndBox,includePartialBlocks)
280 integer, intent(in) :: blockType
281 integer,dimension(MAXBLOCKS),intent(out) :: listOfBlocks
282 integer,intent(out) :: count
283 integer,intent(IN),optional :: refinementLevel
284 real, dimension(LOW:HIGH,MDIM), intent(IN), optional :: region_bndBox
285 logical, intent(IN), optional :: includePartialBlocks
286 end subroutine Grid_getListOfBlocks
287 end interface
288
289 interface
290 subroutine Grid_getLocalNumBlks(numBlocks)
291 integer,intent(out) :: numBlocks
292 end subroutine Grid_getLocalNumBlks
293 end interface
294
295 interface
296 subroutine Grid_getNumBlksFromType(blockType,numBlocks)
297 integer,intent(in) :: blockType
298 integer,intent(out) :: numBlocks
299 end subroutine Grid_getNumBlksFromType
300 end interface
301
303 subroutine Grid_getMinCellSize(minCellSize)
304 real, intent(OUT) :: minCellSize
305 end subroutine Grid_getMinCellSize
306 end interface
307
308 interface
309 subroutine Grid_getMinCellSizes(minCellSizes)
310 real, dimension(MDIM), intent(OUT) :: minCellSizes
311 end subroutine Grid_getMinCellSizes
312 end interface
313
314 interface
315 subroutine Grid_subcellGeometry(nsubI,nsubJ,nsubK, &
316 dvCell, dvSub, xL,xR, yL,yR, pos, blockID)
317 implicit none
318 integer, VALUE_INTENT(IN) :: nsubI, nsubJ, nsubK
319 real, intent(in) :: dvCell
320 real, intent(OUT) :: dvSub(nsubI, nsubJ)
321 real,OPTIONAL,intent(in) :: xL, xR
322 real,OPTIONAL,intent(in) :: yL, yR
323 integer,OPTIONAL, intent(in) :: blockID
324 integer,OPTIONAL, intent(in) :: pos(*)
325 end subroutine Grid_subcellGeometry
326 end interface
327
328
330 subroutine Grid_getSingleCellCoords(ind, level,edge, coords)
331 implicit none
332 integer,dimension(MDIM), intent(in) :: ind
333 integer, intent(in) :: level, edge
334 real, dimension(MDIM), intent(out) :: coords
335 end subroutine Grid_getSingleCellCoords
336 end interface
337
339 subroutine Grid_getSingleCellVol(point, level, cellvolume)
340 integer, intent(in) :: point(1:MDIM)
341 integer, intent(in) :: level
342 real, intent(out) :: cellvolume
343 end subroutine Grid_getSingleCellVol
344 end interface Grid_getSingleCellVol
345
346 interface
347 subroutine Grid_guardCellMaskHook(ccMask, needEos)
348 implicit none
349 logical,intent(INOUT) :: ccMask(*)
350 logical,intent(IN) :: needEos
351 end subroutine Grid_guardCellMaskHook
352 end interface
353
354 interface Grid_init
355 subroutine Grid_init()
356 end subroutine Grid_init
357 end interface
358
360 subroutine Grid_initDomain( restart,particlesInitialized)
361 logical, intent(IN) :: restart
362 logical,intent(INOUT) :: particlesInitialized
363 end subroutine Grid_initDomain
364 end interface
365
366 interface
367 subroutine Grid_makeVector(vecLen,numVars,newVec,numVec,vecLastFree,copyDirection,gridDataStruct)
368 implicit none
369 integer, intent(in) :: vecLen
370 integer, intent(in) :: numVars
371 integer,intent(INOUT) :: numVec
372 real, dimension(vecLen,numVars,numVec),intent(INOUT) :: newVec
373 integer, optional,intent(OUT):: vecLastFree
374 integer, optional,intent(in) :: copyDirection
375 integer, optional,intent(in) :: gridDataStruct
376 end subroutine Grid_makeVector
377 end interface
378
381 end subroutine Grid_markRefineDerefine
382 end interface
383
385 subroutine Grid_markRefineSpecialized(criterion,size,specs,lref)
386 integer, intent(IN) :: criterion
387 integer, intent(IN) :: size
388 real,dimension(size),intent(IN) :: specs
389 integer, intent(IN) :: lref
390 end subroutine Grid_markRefineSpecialized
391 end interface
392
394 subroutine Grid_moveParticles(dataBuf, propCount, maxCount, localCount, &
395 index_list, indexCount,&
396 coords_in_blk)
397
398
399 integer,intent(IN) :: maxCount, propCount, indexCount
400 integer,intent(INOUT) :: localCount
401
402 real, dimension(propCount, maxCount),intent(INOUT) :: dataBuf
403 integer, dimension(indexCount),intent(IN) :: index_list
404 logical, intent(IN) :: coords_in_blk
405
406 end subroutine Grid_moveParticles
407 end interface
408
410 subroutine Grid_notifySolnDataUpdate(gds,mask)
411 implicit none
412 integer,OPTIONAL,intent(in) :: gds !"grid data struct", i.e., CENTER, CENTER_FACES, FACES
413 logical,OPTIONAL,intent(in) :: mask(*) !optional mask, as for Grid_fillGuardCells
414 end subroutine Grid_notifySolnDataUpdate
415 subroutine Grid_notifySolnDataUpdateVlist(varList,gds)
416 implicit none
417 integer,intent(in) :: varList(:) ! list of UNK (or other?) variables
418 integer,OPTIONAL,intent(in) :: gds !"grid data struct", i.e., CENTER, CENTER_FACES, FACES
419 end subroutine Grid_notifySolnDataUpdateVlist
420 end interface
421
422
424 subroutine Grid_putFluxData(level,axis, pressureSlots)
425 implicit none
426 integer, intent(IN) :: level
427 integer, intent(IN),optional :: axis
428 integer, intent(IN), OPTIONAL,target :: pressureSlots(:)
429 end subroutine Grid_putFluxData
430 subroutine Grid_putFluxData_block(blockDesc,fluxBufX,fluxBufY,fluxBufZ, lo,isFluxDensity)
431 use Grid_tile, ONLY : Grid_tile_t
432 implicit none
433 type(Grid_tile_t), intent(in) :: blockDesc
434 integer,intent(in) :: lo(3)
435 real,CONTIGUOUS_INTENT(in),dimension(: ,lo(1): ,lo(2): ,lo(3): ),TARGET :: fluxBufX,fluxBufY,fluxBufZ
436 logical, intent(IN), OPTIONAL :: isFluxDensity(:) !maybe eliminate
437 end subroutine Grid_putFluxData_block
438 end interface Grid_putFluxData
439
440 interface Grid_getFluxData ! get (corrected) flux data
441 subroutine Grid_getFluxData_block(blockDesc,fluxBufX,fluxBufY,fluxBufZ, lo, axis, isFluxDensity)
442 use Grid_tile, ONLY : Grid_tile_t
443 implicit none
444 type(Grid_tile_t), intent(in) :: blockDesc
445 integer,intent(in) :: lo(3)
446 real,intent( OUT),dimension(: ,lo(1): ,lo(2): ,lo(3): ),TARGET :: fluxBufX, fluxBufY, fluxBufZ
447 integer, intent(IN),optional :: axis
448 logical, intent(IN), OPTIONAL :: isFluxDensity(:) !maybe eliminate
449 end subroutine Grid_getFluxData_block
450 end interface
451
453 subroutine Grid_getFluxCorrData_xtra(blockDesc,fluxBufX,fluxBufY,fluxBufZ, lo, fluxCorrX,fluxCorrY,fluxCorrZ, isFluxDensity)
454 use Grid_tile, ONLY : Grid_tile_t
455 implicit none
456 type(Grid_tile_t), intent(in) :: blockDesc
457 integer,intent(in) :: lo(3)
458 real,CONTIGUOUS_INTENT(in) ,dimension(: ,lo(1): ,lo(2): ,lo(3): ),TARGET :: fluxBufX, fluxBufY, fluxBufZ
459 real,CONTIGUOUS_INTENT(OUT),dimension(: ,lo(1): ,lo(2): ,lo(3): ),TARGET :: fluxCorrX,fluxCorrY,fluxCorrZ
460 logical, intent(IN), OPTIONAL :: isFluxDensity(:) !maybe eliminate
461 end subroutine Grid_getFluxCorrData_xtra
462 subroutine Grid_getFluxCorrData_block(blockDesc,fluxBufX,fluxBufY,fluxBufZ, lo, isFluxDensity)
463 use Grid_tile, ONLY : Grid_tile_t
464 implicit none
465 type(Grid_tile_t), intent(in) :: blockDesc
466 integer,intent(in) :: lo(3)
467 real,intent(OUT) ,dimension(: ,lo(1): ,lo(2): ,lo(3): ),TARGET :: fluxBufX, fluxBufY, fluxBufZ
468 logical, intent(IN), OPTIONAL :: isFluxDensity(:) !maybe eliminate
469 end subroutine Grid_getFluxCorrData_block
470 end interface
471
473 subroutine Grid_correctFluxData(blockDesc,fluxBufX,fluxBufY,fluxBufZ, lo, isFluxDensity)
474 use Grid_tile, ONLY : Grid_tile_t
475 implicit none
476 type(Grid_tile_t), intent(in) :: blockDesc
477 integer,intent(in) :: lo(3)
478 real,CONTIGUOUS_INTENT(INOUT),dimension(: ,lo(1): ,lo(2): ,lo(3): ),TARGET :: fluxBufX, fluxBufY, fluxBufZ
479 logical, intent(IN), OPTIONAL :: isFluxDensity(:) !maybe eliminate
480 end subroutine Grid_correctFluxData
481
482 subroutine Grid_correctFluxData_xtra(blockDesc,scaleF,fluxBufX,fluxBufY,fluxBufZ, lo, &
483 scaleC,fluxOldX,fluxOldY,fluxOldZ, &
484 isFluxDensity)
485 use Grid_tile, ONLY : Grid_tile_t
486 implicit none
487 type(Grid_tile_t), intent(in) :: blockDesc
488 integer,intent(in) :: lo(3)
489 real,intent(in) :: scaleF,scaleC
490 real,CONTIGUOUS_INTENT(in ),dimension(: ,lo(1): ,lo(2): ,lo(3): ),TARGET :: fluxOldX, fluxOldY, fluxOldZ
491 real,INTENT(INOUT),dimension(size(fluxOldX,1),size(fluxOldX,2),size(fluxOldX,3),size(fluxOldX,4)),TARGET :: FLUXBUFX !CAPITALIZED TO PREVENT INDEX REORDERING!
492 real,INTENT(INOUT),dimension(size(fluxOldY,1),size(fluxOldY,2),size(fluxOldY,3),size(fluxOldY,4)),TARGET :: FLUXBUFY !CAPITALIZATION INTENTIONAL!
493 real,INTENT(INOUT),dimension(size(fluxOldZ,1),size(fluxOldZ,2),size(fluxOldZ,3),size(fluxOldZ,4)),TARGET :: FLUXBUFZ !CAPITALIZATION INTENTIONAL!
494 logical, intent(IN), OPTIONAL :: isFluxDensity(:) !maybe eliminate
495 end subroutine Grid_correctFluxData_xtra
496 end interface Grid_correctFluxData
497
498
500 subroutine Grid_putLocalNumBlks(numBlocks)
501 integer,intent(in) :: numBlocks
502 end subroutine Grid_putLocalNumBlks
503 end interface
504
507 end subroutine Grid_restrictAllLevels
508 end interface
509
510 interface
511 subroutine Grid_restrictByLevels( gridDataStruct, fromLevel, toLevel, checkFinestLevel,&
512 maskSize,mask)
513 integer, intent(in) :: gridDataStruct
514 integer, intent(in) :: fromLevel, toLevel
515 logical, optional,intent(in) :: checkFinestLevel
516 integer, optional,intent(in) :: maskSize
517 logical,dimension(*),optional,intent(in) :: mask
518 end subroutine Grid_restrictByLevels
519 end interface
520
523 end subroutine Grid_sendOutputData
524 end interface
525
526 interface
527 subroutine Grid_setFluxHandling(handling, status)
528 implicit none
529 character(len=*),intent(IN) :: handling
530 integer,intent(OUT),OPTIONAL :: status
531 end subroutine Grid_setFluxHandling
532 end interface
533
534 interface
535 subroutine Grid_setInterpValsGcell(setval)
536 implicit none
537 logical, intent(IN) :: setval
538 end subroutine Grid_setInterpValsGcell
539 end interface
540
541 interface
542 subroutine Grid_setWork(tileDesc,work,mode)
543 use Grid_tile, ONLY : Grid_tile_t
544 implicit none
545 type(Grid_tile_t),intent(IN) :: tileDesc
546 real, intent(IN) :: work
547 integer, intent(IN),optional :: mode
548 end subroutine Grid_setWork
549 end interface
550
551 interface
552 subroutine Grid_setWorkBounds(pwork_bnd,lwork_bnd,bnd_type)
553 implicit none
554 real, intent(IN) :: pwork_bnd,lwork_bnd
555 integer,intent(IN) :: bnd_type
556 end subroutine Grid_setWorkBounds
557 end interface
558
559 interface
560 subroutine Grid_setWorkDefault(pwork,lwork)
561 implicit none
562 real, intent(IN) :: pwork,lwork
563 end subroutine Grid_setWorkDefault
564 end interface
565
567 subroutine Grid_updateRefinement( nstep,time, gridChanged)
568 integer, intent(in) :: nstep
569 real, intent(in) :: time
570 logical, intent(out), OPTIONAL :: gridChanged
571 end subroutine Grid_updateRefinement
572 end interface
573
575 subroutine Grid_unitTest(fileUnit,perfect)
576
577 integer, intent(in) :: fileUnit ! Output to file
578 logical, intent(inout) :: perfect ! Flag to indicate errors
579 end subroutine Grid_unitTest
580 end interface
581
583 subroutine Grid_getGeometry(geometry)
584 integer, intent(OUT) :: geometry
585 end subroutine Grid_getGeometry
586 end interface
587
588
590 subroutine Grid_sortParticles(dataBuf,props,localNumCount,&
591 elementTypes,maxPerProc,&
592 elementsPerBlk, attrib1, attrib2)
593
594 implicit none
595 integer,intent(INOUT) :: localNumCount
596 integer,intent(IN) :: maxPerProc, props,elementTypes
597 real,intent(INOUT),dimension(props,maxPerProc) :: dataBuf
598 integer,intent(OUT),dimension(MAXBLOCKS,elementTypes) :: elementsPerBlk
599 integer, intent(IN) :: attrib1
600 integer, optional, intent(IN) :: attrib2
601 end subroutine Grid_sortParticles
602 end interface
603
605 subroutine Grid_countParticles(props, localCount,elementTypes,maxCount)
606
607 implicit none
608 integer,intent(INOUT) :: localCount
609 integer,intent(IN) :: elementTypes
610 integer, optional, intent(IN) :: maxCount
611 integer,intent(IN) :: props
612 end subroutine Grid_countParticles
613 end interface
614
616 subroutine Grid_countParticlesByBlock(particlesPerBlock)
617
618 implicit none
619 integer,dimension(MAXBLOCKS, NPART_TYPES), intent(OUT) :: particlesPerBlock
620 end subroutine Grid_countParticlesByBlock
621 end interface
622
623
625 subroutine Grid_mapMeshToParticles_pc (ptContainerPos, part_props,part_blkID,&
626 posAttrib,&
627 numAttrib, attrib,&
628 mapType,gridDataStruct)
629 implicit none
630 integer, INTENT(in) :: ptContainerPos, part_props, part_blkID
631 integer, intent(IN) :: numAttrib
632 integer, dimension(PART_ATTR_DS_SIZE,numAttrib),INTENT(in) :: attrib
633 integer,dimension(MDIM), intent(IN) :: posAttrib
634 integer, INTENT(IN) :: mapType
635 integer, optional, intent(IN) :: gridDataStruct
636 end subroutine Grid_mapMeshToParticles_pc
637
638 subroutine Grid_mapMeshToParticles (particles, part_props,part_blkID,&
639 numParticles,&
640 posAttrib,&
641 numAttrib, attrib,&
642 mapType,gridDataStruct)
643 implicit none
644 integer, INTENT(in) :: part_props, numParticles, part_blkID
645 real, INTENT(inout),dimension(part_props,numParticles) :: particles
646 integer, intent(IN) :: numAttrib
647 integer, dimension(PART_ATTR_DS_SIZE,numAttrib),INTENT(in) :: attrib
648 integer,dimension(MDIM), intent(IN) :: posAttrib
649 integer, INTENT(IN) :: mapType
650 integer, optional, intent(IN) :: gridDataStruct
651 end subroutine Grid_mapMeshToParticles
652 end interface
653
654 interface
655 subroutine Grid_mapParticlesToMesh (particles,part_props,numParticles,&
656 maxParticlesPerProc,propPart, varGrid, mode, ptInfo)
657
658 integer,intent(IN) :: numParticles, part_props,maxParticlesPerProc
659 real,dimension(part_props,numParticles),intent(INOUT) :: particles
660 integer, INTENT(in) :: propPart, varGrid
661 integer, INTENT(in), optional :: mode
662 integer, INTENT(in), optional :: ptInfo
663 end subroutine Grid_mapParticlesToMesh
664 end interface
665
666
667 interface
668 subroutine Grid_solvePoisson (iSoln, iSrc, bcTypes, &
669 bcValues, poisfact)
670 implicit none
671 integer, intent(in) :: iSoln, iSrc
672 integer, intent(in) :: bcTypes(6)
673 real, intent(in) :: bcValues(2,6)
674 real, intent(inout) :: poisfact
675 end subroutine Grid_solvePoisson
676 end interface
677
678 interface
679 subroutine Grid_setSolverDbgContextInfo(component,group)
680 implicit none
681 integer,intent(in),OPTIONAL :: component, group
682 end subroutine Grid_setSolverDbgContextInfo
683 end interface
684
685 interface
686 subroutine Grid_limitAbundance(blkLimits,solnData)
687 integer, dimension(2,MDIM), INTENT(in) :: blkLimits
688 real, POINTER :: solnData(:,:,:,:)
689 end subroutine Grid_limitAbundance
690 end interface
691
692
693 interface
694 subroutine Grid_renormAbundance(tileDesc, tileLimits, solnData)
695 use Grid_tile, ONLY : Grid_tile_t
696 implicit none
697 type(Grid_tile_t), intent(IN) :: tileDesc
698 integer, intent(IN) :: tileLimits(LOW:HIGH, 1:MDIM)
699 real, pointer :: solnData(:,:,:,:)
700 end subroutine Grid_renormAbundance
701 end interface
702
703 interface
704 subroutine Grid_renormMassScalars(blkLimits,solnData)
705 integer, intent(in), dimension(2,MDIM)::blkLimits
706 real,pointer :: solnData(:,:,:,:)
707 end subroutine Grid_renormMassScalars
708 end interface
709
710 interface
711 subroutine Grid_conserveField ()
712 end subroutine Grid_conserveField
713 end interface
714
715
716 interface
717 subroutine Grid_pfft(direction,inArray,outArray)
718 integer, intent(IN) :: direction
719 real, dimension(:),intent(IN) :: inArray
720 real,dimension(:), intent(OUT) :: outArray
721 end subroutine Grid_pfft
722 end interface
723
724 interface
725 subroutine Grid_pfftInit( ndim, needMap,globalLen, mapSize,&
726 transformType, baseDatType, jProcs, kProcs, refinementLevel, region_bndBox)
727 integer, intent(IN) :: ndim
728 logical, intent(IN) :: needMap
729 integer, dimension(MDIM), intent(IN) :: globalLen
730 integer,dimension(MDIM),intent(OUT) :: mapSize
731 integer,dimension(MDIM),optional,intent(IN) :: transformType
732 integer,dimension(0:MDIM),optional,intent(IN) :: baseDatType
733 integer,optional,intent(IN) :: jProcs, kProcs
734 integer,optional,intent(IN) :: refinementLevel
735 real, dimension(LOW:HIGH,MDIM), optional, intent(IN) :: region_bndBox
736 end subroutine Grid_pfftInit
737 end interface
738
739 interface
740 subroutine Grid_pfftFinalize()
741 end subroutine Grid_pfftFinalize
742 end interface
743
744 interface
745 subroutine Grid_bcApplyToRegionSpecialized(bcType,gridDataStruct,level,&
746 guard,axis,face,regionData,regionSize,mask,applied,&
747 secondDir,ThirdDir,endPoints,idest)
748 implicit none
749
750 integer, intent(IN) :: bcType,axis,face,guard,gridDataStruct, level
751 integer,dimension(REGION_DIM),intent(IN) :: regionSize
752 real,dimension(regionSize(BC_DIR),&
753 regionSize(SECOND_DIR),&
754 regionSize(THIRD_DIR),&
755 regionSize(STRUCTSIZE)),intent(INOUT)::regionData
756 logical,intent(IN),dimension(regionSize(STRUCTSIZE)):: mask
757 logical, intent(OUT) :: applied
758 integer,intent(IN) :: secondDir,thirdDir
759 integer,intent(IN),dimension(LOW:HIGH,MDIM) :: endPoints
760 integer,intent(IN),OPTIONAL:: idest
762 end interface
763
764 interface
765 subroutine Grid_bcApplyToRegion(bcType,gridDataStruct,level,&
766 guard,axis,face,regionData,regionSize,mask,applied,&
767 secondDir,ThirdDir,endPoints,idest)
768 implicit none
769 integer, intent(IN) :: bcType,axis,face,guard,gridDataStruct, level
770 integer,dimension(REGION_DIM),intent(IN) :: regionSize
771 real,dimension(regionSize(BC_DIR),&
772 regionSize(SECOND_DIR),&
773 regionSize(THIRD_DIR),&
774 regionSize(STRUCTSIZE)),intent(INOUT)::regionData
775 logical,intent(IN),dimension(regionSize(STRUCTSIZE)):: mask
776 logical, intent(OUT) :: applied
777 integer,intent(IN) :: secondDir,thirdDir
778 integer,intent(IN),dimension(LOW:HIGH,MDIM) :: endPoints
779 integer,intent(IN),OPTIONAL:: idest
780 end subroutine Grid_bcApplyToRegion
781 end interface
782
783 interface
784 subroutine Grid_bcApplyToRegionMixedGds(bcType,gridDataStruct,level,&
785 guard,axis,face,&
786 regionDataC,regionDataFN,regionDataFT1,regionDataFT2,&
787 regionSizeCtr,&
788 applied,&
789 tileDesc,secondDir,thirdDir,endPointsCtr,rightHanded,idest)
790 use Grid_tile, ONLY : Grid_tile_t
791 implicit none
792 integer, intent(IN) :: bcType,axis,face,guard,gridDataStruct,level
793 integer,dimension(REGION_DIM),intent(IN) :: regionSizeCtr
794 real,pointer,dimension(:,:,:,:) :: regionDataFN, regionDataFT1, regionDataFT2, regionDataC
795 logical, intent(INOUT) :: applied
796 type(Grid_tile_t),intent(IN) :: tileDesc
797 integer,intent(IN) :: secondDir,thirdDir
798 integer,intent(IN),dimension(LOW:HIGH,MDIM) :: endPointsCtr
799 logical, intent(IN) :: rightHanded
800 integer,intent(IN),OPTIONAL:: idest
801 end subroutine Grid_bcApplyToRegionMixedGds
802 end interface
803
804 interface
805 subroutine Grid_pfftGetIndexLimits(configLimits,phaseLimits)
806 integer,dimension(LOW:HIGH,MDIM),intent(OUT) :: configLimits, phaseLimits
807 end subroutine Grid_pfftGetIndexLimits
808 end interface
809
810
811 interface
812 subroutine Grid_getBlkType(blockID, blkType)
813 implicit none
814 integer,intent(in) :: blockID
815 integer,intent(out) :: blkType
816 end subroutine Grid_getBlkType
817 end interface
818
820 subroutine Grid_pfftMapToInput(gridVar, pfft_inArray)
821 integer,intent(IN) :: gridVar
822 real, dimension(:),intent(OUT) :: pfft_inArray
823 end subroutine Grid_pfftMapToInput
824 subroutine Grid_pfftMapToInput3DArr(gridVar, pfft_inArray)
825 integer,intent(IN) :: gridVar
826 real, dimension(:,:,:),intent(OUT) :: pfft_inArray
827 end subroutine Grid_pfftMapToInput3DArr
828 end interface
829
831 subroutine Grid_pfftMapFromOutput(gridVar, pfft_outArray)
832 integer,intent(IN) :: gridVar
833 real, dimension(:),intent(IN) :: pfft_outArray
834 end subroutine Grid_pfftMapFromOutput
835 subroutine Grid_pfftMapFromOutput3DArr(gridVar, pfft_outArray)
836 integer,intent(IN) :: gridVar
837 real, dimension(:,:,:),intent(IN) :: pfft_outArray
838 end subroutine Grid_pfftMapFromOutput3DArr
839 end interface
840
841 interface
842 subroutine Grid_outsideBoundBox(pos,bndBox,outside,Negh)
843 real,dimension(MDIM),intent(IN) :: pos
844 real,dimension(LOW:HIGH,MDIM),intent(IN) :: bndBox
845 logical, intent(OUT) :: outside
846 integer, dimension(MDIM),intent(OUT) :: Negh
847 end subroutine Grid_outsideBoundBox
848 end interface
849
850 interface
851 subroutine Grid_getMaxCommonRefinement(inputComm, maxRefinement)
852 implicit none
853 integer, intent(IN) :: inputComm
854 integer, intent(OUT) :: maxRefinement
855 end subroutine Grid_getMaxCommonRefinement
856 end interface
857
858 interface
859 subroutine Grid_getMaxRefinement(maxRefinement, mode, scope, inputComm)
860 implicit none
861 integer, intent(IN), OPTIONAL :: mode, scope
862 integer, intent(IN), OPTIONAL :: inputComm
863 integer, intent(OUT) :: maxRefinement
864 end subroutine Grid_getMaxRefinement
865 end interface
866
867 interface
868
869 subroutine Grid_GCPutScratch(gridDataStruct,needSetup,releaseSetup,&
870 &blkCount,blkList,indCount,indList, gcCnt, putData)
871
872 integer, intent(IN) :: gridDataStruct
873 logical, intent(IN) :: needSetup, releaseSetup
874 integer, optional,intent(IN) :: blkCount, indCount
875 integer, optional, dimension(:), intent(IN) :: blkList
876 integer, optional, dimension(:), intent(IN) :: indList
877 integer, optional, dimension(NDIM),intent(IN) :: gcCnt
878 logical, optional,intent(IN) :: putData
879
880 end subroutine Grid_GCPutScratch
881
882 end interface
883
884 interface
885
886 subroutine Grid_GCTransferOneBlk(mode,gridDataStruct,blkIndex,blkData)
887
888 logical, intent(IN) :: mode
889 integer, intent(IN) :: gridDataStruct, blkIndex
890 real, pointer, dimension(:,:,:,:) :: blkData
891
892 end subroutine Grid_GCTransferOneBlk
893
894 end interface
895
896
897 interface
898 subroutine Grid_getNumVars(gridStruct, nVar)
899 implicit none
900 integer, intent(in) :: gridStruct
901 integer, intent(out) :: nVar
902 end subroutine Grid_getNumVars
903 end interface
904
905 interface
906
907 subroutine Grid_addToVar(srcVar, destVar, multFactor, reset)
908 integer, intent(in) :: srcVar, destVar
909 real, intent(in) :: multFactor
910 logical, intent(in) :: reset
911
912 end subroutine Grid_addToVar
913
914 end interface
915
916 interface
917 subroutine Grid_smoothVar(ivar, ivarOut, &
918 solnData, lbUI,lbUJ,lbUK, &
919 blklim, smoothMethod, gcLayers, blockID,&
920 useMinSmoothVarVal,&
921 minSmoothVarVal,&
922 useMaxSmoothVarVal,&
923 maxSmoothVarVal,&
924 smoothCoeff )
925 integer, intent(in) :: ivar
926 integer, intent(in) :: ivarOut
927 integer, VALUE_INTENT(IN) :: lbUI,lbUJ,lbUK
928 real, intent(INOUT) :: solnData(:,lbUI:,lbUJ:,lbUK:)
929 integer, intent(in) :: blklim(2,MDIM)
930 integer, intent(in),OPTIONAL :: smoothMethod
931 integer, intent(IN),OPTIONAL :: gcLayers
932 integer, intent(IN),OPTIONAL :: blockID
933 logical, intent(IN),OPTIONAL :: useMinSmoothVarVal,useMaxSmoothVarVal
934 real , intent(IN),OPTIONAL :: minSmoothVarVal,maxSmoothVarVal
935 real , intent(IN),OPTIONAL :: smoothCoeff
936 end subroutine Grid_smoothVar
937 end interface
938
939 interface
940 subroutine Grid_primitiveToConserve(blkList,count,force)
941 integer,intent(IN) :: count
942 integer,dimension(count),intent(IN) :: blkList
943 logical,intent(IN) :: force
944 end subroutine Grid_primitiveToConserve
945 end interface
946
947 interface
948 subroutine Grid_conserveToPrimitive(blkList,count,allCells,force)
949 integer,intent(IN) :: count
950 integer,dimension(count),intent(IN) :: blkList
951 logical,intent(IN) :: allCells,force
952 end subroutine Grid_conserveToPrimitive
953 end interface
954
955 interface
956
957 subroutine Grid_getDomainBoundBox(boundBox)
958 real,dimension(LOW:HIGH,MDIM), intent(OUT) :: boundBox
959 end subroutine Grid_getDomainBoundBox
960 end interface
961
962 interface
963 subroutine Grid_getDomainBC(boundary)
964 integer,dimension(LOW:HIGH,MDIM), intent(OUT) :: boundary
965 end subroutine Grid_getDomainBC
966 end interface
967
968 interface
969 subroutine Grid_parseNonRep(strlwr, nonrep, idx)
970 implicit none
971 character(*), intent(in) :: strlwr
972 integer, intent(out) :: nonrep, idx
973 end subroutine
974 end interface
975
976 interface
977 subroutine Grid_formatNonRep(nonrep, idx, str)
978 implicit none
979 integer, intent(in) :: nonrep, idx
980 character(*), intent(out) :: str
981 end subroutine
982 end interface
983
984 interface
985 subroutine Grid_getVarNonRep(mapblock, var, nonrep, locidx)
986 implicit none
987 integer, intent(in) :: mapblock, var
988 integer, intent(out) :: nonrep
989 integer, intent(out), optional :: locidx
990 end subroutine
991 end interface
992
994
995 subroutine Grid_getBlkNeighBlkIDFromPos(blockDesc, pos, neghDir, ansBlockID, ansProcID)
996 use Grid_tile, ONLY : Grid_tile_t
997 implicit none
998 type(Grid_tile_t), intent(IN) :: blockDesc
999 real, dimension(MDIM), intent(IN) :: pos
1000 integer, dimension(MDIM), intent(IN) :: neghDir
1001 integer, intent(OUT) :: ansBlockID, ansProcID
1002 end subroutine Grid_getBlkNeighBlkIDFromPos
1003
1004 subroutine Grid_getBlkIDFromPos(pos, ansBlockID, ansProcID, comm)
1005 implicit none
1006 real, dimension(1:MDIM), intent(IN) :: pos
1007 integer, intent(OUT) :: ansBlockID, ansProcID
1008 integer,optional,intent(IN) :: comm
1009 end subroutine Grid_getBlkIDFromPos
1010 end interface
1011
1013 subroutine Grid_getLocalBlkIDFromPosSimple(pos,ansBlockID, ansProcID,blkList, blkCount,blockType)
1014 implicit none
1015 real, dimension(1:MDIM), intent(IN) :: pos
1016 integer, intent(OUT) :: ansBlockID, ansProcID
1017 integer, OPTIONAL,intent(IN) :: blkCount
1018 integer, OPTIONAL,dimension(:),intent(IN),target :: blkList
1019 integer, OPTIONAL, intent(IN) :: blockType
1020 end subroutine Grid_getLocalBlkIDFromPosSimple
1021 end interface
1022
1023 interface
1025 end subroutine Grid_sbBroadcastParticles
1026 end interface
1027
1028 interface
1030 end subroutine Grid_sbCreateGroups
1031 end interface
1032
1033 interface
1035 end subroutine Grid_sbSelectMaster
1036 end interface
1037
1038! interface
1039! subroutine Grid_updateSolidBodyForces()
1040! end subroutine Grid_updateSolidBodyForces
1041! end interface
1042
1043 interface
1044 subroutine Grid_updateSolidBodyForces(blkID,particleData)
1045 implicit none
1046 integer, intent(in) :: blkID
1047 real, intent(inout) :: particleData(NPART_PROPS)
1048 end subroutine Grid_updateSolidBodyForces
1049 end interface
1050
1051
1052 interface
1053 subroutine Grid_solidBodyUnitTest(fileUnit, perfect)
1054 integer, intent(IN) :: fileUnit
1055 logical, intent(INOUT) :: perfect
1056 end subroutine Grid_solidBodyUnitTest
1057 end interface
1058
1059 interface
1060 subroutine Grid_receiveInputData(localNumBlocks, alnblocks, xx)
1061 implicit none
1062 integer, intent(IN) :: localNumBlocks, alnblocks, xx
1063 end subroutine Grid_receiveInputData
1064 end interface
1065
1066 interface
1067 subroutine Grid_getNeighProcList(includeMyProc, neighProcList, numNeigh)
1068 implicit none
1069 logical, intent(IN) :: includeMyProc
1070 integer, dimension(:), pointer :: neighProcList
1071 integer, intent(OUT) :: numNeigh
1072 end subroutine Grid_getNeighProcList
1073 end interface
1074
1075 interface
1076 subroutine Grid_getBlkNeighLevels(blockID, levels, trackBdry)
1077 implicit none
1078 integer,intent(in) :: blockID
1079 integer,intent(OUT) :: levels(-1:1, -K2D:K2D , -K3D:K3D)
1080 logical,intent(in),OPTIONAL :: trackBdry
1081 end subroutine Grid_getBlkNeighLevels
1082 end interface
1083
1084 interface
1085 subroutine Grid_coordTransfm(x,y,z, xout,yout,zout, geometryIn,geometryOut, ndim, velI,velJ,velK,velIOut,velJOut,velKOut)
1086 implicit none
1087 real,intent(IN) :: x,y,z
1088 real,intent(OUT) :: xout,yout,zout
1089 integer,OPTIONAL,intent(IN) :: geometryIn, geometryOut
1090 integer,OPTIONAL,intent(IN) :: ndim
1091 real,OPTIONAL,intent(IN) :: velI,velJ,velK
1092 real,OPTIONAL,intent(OUT) :: velIOut,velJOut,velKOut
1093 end subroutine Grid_coordTransfm
1094 end interface
1095
1096
1097 interface
1098 subroutine Grid_getTileIterator(itor, nodetype, level, tiling, tileSize, nthreads)
1099 use Grid_iterator, ONLY : Grid_iterator_t
1100 implicit none
1101 type(Grid_iterator_t), intent(OUT) :: itor
1102 integer, intent(IN) :: nodetype
1103 integer, intent(IN), optional :: level
1104 logical, intent(IN), optional :: tiling
1105 integer, intent(IN), optional :: tileSize(1:MDIM)
1106 integer, intent(IN), optional :: nthreads
1107 end subroutine Grid_getTileIterator
1108 end interface
1109
1110 interface
1112 use Grid_iterator, ONLY : Grid_iterator_t
1113 implicit none
1114 type(Grid_iterator_t), intent(INOUT) :: itor
1115 end subroutine Grid_releaseTileIterator
1116 end interface
1117
1118 interface
1119 logical function Grid_properTilingWanted()
1120 implicit none
1121 end function Grid_properTilingWanted
1122 end interface
1123
1124 interface
1126 implicit none
1127 end subroutine Grid_zeroFluxData
1128 end interface
1129
1130 interface
1131 subroutine Grid_addFineToFluxRegister(level, isDensity, coefficient, &
1132 zeroFullRegister)
1133 implicit none
1134 integer, intent(IN) :: level
1135 logical, intent(IN), optional :: isDensity(:)
1136 real, intent(IN), optional :: coefficient
1137 logical, intent(IN), optional :: zeroFullRegister
1138 end subroutine Grid_addFineToFluxRegister
1139 end interface
1140
1141 interface
1142 subroutine Grid_addCoarseToFluxRegister(level, isDensity, coefficient, &
1143 zeroFullRegister)
1144 implicit none
1145 integer, intent(IN) :: level
1146 logical, intent(IN), optional :: isDensity(:)
1147 real, intent(IN), optional :: coefficient
1148 logical, intent(IN), optional :: zeroFullRegister
1149 end subroutine Grid_addCoarseToFluxRegister
1150 end interface
1151
1152end Module Grid_interface
1153
subroutine Grid_ascGetBlk5Ptr(blockID, data5Ptr, gridDataStruct)
subroutine Grid_correctFluxData_xtra(blockDesc, scaleF, fluxBufX, fluxBufY, fluxBufZ, lo, scaleC, fluxOldX, fluxOldY, fluxOldZ, isFluxDensity)
subroutine Grid_getBlkNeighBlkIDFromPos(blockDesc, pos, neghDir, ansBlockID, ansProcID)
subroutine Grid_getFluxCorrData_block(blockDesc, fluxBufX, fluxBufY, fluxBufZ, lo, isFluxDensity)
subroutine Grid_getFluxCorrData_xtra(blockDesc, fluxBufX, fluxBufY, fluxBufZ, lo, fluxCorrX, fluxCorrY, fluxCorrZ, isFluxDensity)
subroutine Grid_getFluxData_block(blockDesc, fluxBufX, fluxBufY, fluxBufZ, lo, axis, isFluxDensity)
subroutine Grid_getLocalBlkIDFromPosSimple(pos, ansBlockID, ansProcID, blkList, blkCount, blockType)
subroutine Grid_mapMeshToParticles_pc(pt_containerPos, part_props, part_blkID, posAttrib, numAttrib, attrib, mapType, gridDataStruct)
subroutine Grid_notifySolnDataUpdateVlist(varList, gds)
subroutine Grid_putFluxData_block(blockDesc, fluxBufX, fluxBufY, fluxBufZ, lo, isFluxDensity)
subroutine Grid_ascDeallocMem(gds, arrayRank)
subroutine Grid_ascStart()
subroutine Grid_ascAllocMem(gds, var1, nvars, nGuardCtr, nGuardFaceN, nGuardFaceT, leafBlocks, arrayRank, highSize)
subroutine Grid_getBlkIndexLimits(blockId, blkLimits, blkLimitsGC, gridDataStruct)
integer, parameter GRID_PDE_BND_DIRICHLET
integer, parameter GRID_PDE_BND_ISOLATED
integer, parameter GRID_PDE_BND_GIVENGRAD
integer, parameter GRID_COPYDIR_TO_VECT
integer, parameter GRID_PDE_BND_NEUMANN
integer, parameter GRID_PDE_BND_PERIODIC
integer, parameter GRID_PDE_BND_GIVENVAL
integer, parameter GRID_COPYDIR_FROM_VECT