grid_structured.F90       coverage:  71.43 %func     59.51 %block


     1) module Grid_Structured_module
     2) 
     3)   use PFLOTRAN_Constants_module
     4) 
     5)   implicit none
     6)  
     7)   private
     8)  
     9) #include "petsc/finclude/petscsys.h"
    10) 
    11) ! structured grid faces
    12)   PetscInt, parameter, public :: NULL_FACE = 0
    13)   PetscInt, parameter, public :: WEST_FACE = 1
    14)   PetscInt, parameter, public :: EAST_FACE = 2
    15)   PetscInt, parameter, public :: SOUTH_FACE = 3
    16)   PetscInt, parameter, public :: NORTH_FACE = 4
    17)   PetscInt, parameter, public :: BOTTOM_FACE = 5
    18)   PetscInt, parameter, public :: TOP_FACE = 6
    19) 
    20)   PetscInt, parameter, public :: CARTESIAN_GRID = 3
    21)   PetscInt, parameter, public :: CYLINDRICAL_GRID = 4
    22)   PetscInt, parameter, public :: SPHERICAL_GRID = 5
    23) 
    24)   type, public :: grid_structured_type
    25) 
    26)     character(len=MAXWORDLENGTH) :: ctype
    27)     PetscInt :: itype  ! type of grid (e.g. structured, unstructured, etc.)
    28)     PetscInt :: global_offset
    29)     PetscInt :: nx, ny, nz    ! Global domain dimensions of the grid.
    30)     PetscInt :: nxy, nmax     ! nx * ny, nx * ny * nz
    31)     PetscInt :: npx, npy, npz ! Processor partition in each direction.
    32)     PetscInt :: npx_final, npy_final, npz_final ! actual decomposition
    33)     PetscInt :: nlx, nly, nlz ! Local grid dimension w/o ghost nodes.
    34)     PetscInt :: ngx, ngy, ngz ! Local grid dimension with ghost nodes.
    35)     PetscInt :: lxs, lys, lzs 
    36)       ! Global indices of non-ghosted corner (starting) of local domain.
    37)     PetscInt :: gxs, gys, gzs
    38)       ! Global indices of ghosted starting corner of local domain.
    39)     PetscInt :: lxe, lye, lze, gxe, gye, gze
    40)       ! Global indices of non-ghosted/ghosted ending corner of local domain.
    41)     PetscInt :: nlxy, nlxz, nlyz
    42)     PetscInt :: ngxy, ngxz, ngyz
    43)     
    44)     PetscInt :: istart, jstart, kstart, iend, jend, kend
    45)       ! istart gives the ghosted local x-index of the non-ghosted starting 
    46)       ! (lower left)corner. iend gives the local x-index of the non-ghosted 
    47)       ! ending corner. jstart, jend correspond to y-index, kstart, kend to 
    48)       ! z-index.  These are all zero-based indexing    
    49) 
    50)     PetscInt :: nlmax  ! Total number of non-ghosted nodes in local domain.
    51)     PetscInt :: ngmax  ! Number of ghosted & non-ghosted nodes in local domain.
    52)     PetscInt :: nlmax_faces  ! Total number of non-ghosted faces in local domain.
    53)     PetscInt :: ngmax_faces  ! Number of ghosted & non-ghosted faces in local domain.
    54) 
    55)     PetscReal :: local_origin(3) ! local origin of non-ghosted grid
    56)     PetscReal :: bounds(3,2)
    57) 
    58)     ! grid spacing for each direction for global domain
    59)     PetscReal, pointer :: dx_global(:), dy_global(:), dz_global(:)
    60)     ! grid spacing for each direction for local, ghosted domain
    61)     PetscReal, pointer :: dxg_local(:), dyg_local(:), dzg_local(:)
    62)     
    63)     PetscBool :: invert_z_axis
    64)     
    65)     PetscReal, pointer :: dx(:), dy(:), dz(:)  ! ghosted grid spacings for each grid cell
    66)     
    67)     PetscInt, pointer :: cell_neighbors(:,:)
    68)     
    69)   end type grid_structured_type
    70) 
    71)   public :: StructGridCreate, &
    72)             StructGridDestroy, &
    73)             StructGridCreateDM, &
    74)             StructGridComputeLocalBounds, &
    75)             StructGridComputeInternConnect, &
    76)             StructGridCreateVecFromDM, &
    77)             StructGridMapIndices, &
    78)             StructGridComputeSpacing, &
    79)             StructGridComputeCoord, &
    80)             StructGridReadDXYZ, &
    81)             StructGridComputeVolumes, &
    82)             StructGridPopulateConnection, &
    83)             StructGridGetIJKFromCoordinate, &
    84)             StructGridGetIJKFromLocalID, &
    85)             StructGridGetIJKFromGhostedID, &
    86)             StructGridGetGhostedNeighbors, &
    87)             StructGridCreateTVDGhosts, &
    88)             StructGridGetGhostedNeighborsCorners
    89)             
    90) contains
    91) 
    92) ! ************************************************************************** !
    93) 
    94) function StructGridCreate()
    95)   ! 
    96)   ! Creates a structured grid object
    97)   ! 
    98)   ! Author: Glenn Hammond
    99)   ! Date: 10/22/07
   100)   ! 
   101) 
   102)   implicit none
   103)   
   104)   type(grid_structured_type), pointer :: StructGridCreate
   105) 
   106)   type(grid_structured_type), pointer :: structured_grid
   107) 
   108)   allocate(structured_grid)
   109)   
   110)   structured_grid%ctype = ''
   111)   structured_grid%itype = 0
   112)   structured_grid%global_offset = 0
   113)   structured_grid%nx = 0
   114)   structured_grid%ny = 0
   115)   structured_grid%nz = 0
   116)   structured_grid%npx = PETSC_DECIDE
   117)   structured_grid%npy = PETSC_DECIDE
   118)   structured_grid%npz = PETSC_DECIDE
   119)   
   120)   structured_grid%npx_final = 0
   121)   structured_grid%npy_final = 0
   122)   structured_grid%npz_final = 0
   123) 
   124)   structured_grid%nx = 0
   125)   structured_grid%ny = 0
   126)   structured_grid%nz = 0
   127)   
   128)   structured_grid%nxy = 0
   129)   structured_grid%nmax = 0
   130)   
   131)   structured_grid%nlx = 0
   132)   structured_grid%nly = 0
   133)   structured_grid%nlz = 0
   134)   structured_grid%nlxy = 0
   135)   structured_grid%nlxz = 0
   136)   structured_grid%nlyz = 0
   137)   structured_grid%nlmax = 0
   138)   
   139)   structured_grid%ngx = 0
   140)   structured_grid%ngy = 0
   141)   structured_grid%ngz = 0
   142)   structured_grid%ngxy = 0
   143)   structured_grid%ngxz = 0
   144)   structured_grid%ngyz = 0
   145)   structured_grid%ngmax = 0
   146) 
   147)   structured_grid%lxs = 0
   148)   structured_grid%lys = 0
   149)   structured_grid%lzs = 0
   150) 
   151)   structured_grid%gxs = 0
   152)   structured_grid%gys = 0
   153)   structured_grid%gzs = 0
   154) 
   155)   structured_grid%lxe = 0
   156)   structured_grid%lye = 0
   157)   structured_grid%lze = 0
   158) 
   159)   structured_grid%gxe = 0
   160)   structured_grid%gye = 0
   161)   structured_grid%gze = 0
   162) 
   163)   structured_grid%istart = 0
   164)   structured_grid%jstart = 0
   165)   structured_grid%kstart = 0
   166)   structured_grid%iend = 0
   167)   structured_grid%jend = 0
   168)   structured_grid%kend = 0
   169)   
   170)   nullify(structured_grid%dx_global)
   171)   nullify(structured_grid%dy_global)
   172)   nullify(structured_grid%dz_global)
   173) 
   174)   nullify(structured_grid%dxg_local)
   175)   nullify(structured_grid%dyg_local)
   176)   nullify(structured_grid%dzg_local)
   177) 
   178)   nullify(structured_grid%dx)
   179)   nullify(structured_grid%dy)
   180)   nullify(structured_grid%dz)
   181)   
   182)   nullify(structured_grid%cell_neighbors)
   183)  
   184)   structured_grid%local_origin = -1.d20
   185)   structured_grid%bounds = -1.d20
   186)   
   187)   structured_grid%invert_z_axis = PETSC_FALSE
   188)   
   189)   StructGridCreate => structured_grid
   190)   
   191) end function StructGridCreate
   192) 
   193) ! ************************************************************************** !
   194) 
   195) subroutine StructGridCreateDM(structured_grid,da,ndof,stencil_width, &
   196)                               stencil_type,option)
   197)   ! 
   198)   ! StructGridCreateDMs: Creates structured distributed, parallel meshes/grids
   199)   ! 
   200)   ! Author: Glenn Hammond
   201)   ! Date: 10/22/07
   202)   ! 
   203) 
   204)   use Option_module
   205)         
   206)   implicit none
   207) 
   208) #include "petsc/finclude/petscvec.h"
   209) #include "petsc/finclude/petscvec.h90"
   210) #include "petsc/finclude/petscdm.h"
   211) #include "petsc/finclude/petscdm.h90"
   212) #include "petsc/finclude/petscdmda.h"
   213) 
   214)   type(option_type) :: option
   215)   type(grid_structured_type) :: structured_grid
   216)   DM :: da
   217)   PetscInt :: ndof
   218)   PetscInt :: stencil_width
   219)   PetscEnum :: stencil_type
   220) 
   221)   PetscErrorCode :: ierr
   222) 
   223)   !-----------------------------------------------------------------------
   224)   ! Generate the DM object that will manage communication.
   225)   !-----------------------------------------------------------------------
   226)   ! This code is for the DMDACreate3D() interface in PETSc versions >= 3.2 --RTM
   227)   call DMDACreate3D(option%mycomm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE, &
   228)                     DM_BOUNDARY_NONE,stencil_type, &
   229)                     structured_grid%nx,structured_grid%ny,structured_grid%nz, &
   230)                     structured_grid%npx,structured_grid%npy, &
   231)                     structured_grid%npz,ndof,stencil_width, &
   232)                     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
   233)                     da,ierr);CHKERRQ(ierr)
   234)   call DMDAGetInfo(da,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
   235)                    PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
   236)                    structured_grid%npx_final,structured_grid%npy_final, &
   237)                    structured_grid%npz_final, &
   238)                    PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
   239)                    PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
   240)                    PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr)
   241) 
   242) end subroutine StructGridCreateDM
   243) 
   244) ! ************************************************************************** !
   245) 
   246) subroutine StructGridComputeLocalBounds(structured_grid,da,option)
   247)   ! 
   248)   ! Computes the corners for the local portion
   249)   ! of the structured grid
   250)   ! 
   251)   ! Author: Glenn Hammond
   252)   ! Date: 03/13/08
   253)   ! 
   254) 
   255)   use Option_module
   256)   
   257)   implicit none
   258)      
   259) #include "petsc/finclude/petscvec.h"
   260) #include "petsc/finclude/petscvec.h90"
   261) #include "petsc/finclude/petscdm.h"
   262) #include "petsc/finclude/petscdm.h90"
   263) 
   264)   type(grid_structured_type) :: structured_grid
   265)   type(option_type) :: option
   266)   DM :: da
   267) 
   268)   PetscErrorCode :: ierr
   269) 
   270)   ! get corner information
   271)   call DMDAGetCorners(da, structured_grid%lxs, &
   272)       structured_grid%lys, structured_grid%lzs, structured_grid%nlx, &
   273)       structured_grid%nly, structured_grid%nlz, ierr);CHKERRQ(ierr)
   274)      
   275)   structured_grid%lxe = structured_grid%lxs + structured_grid%nlx
   276)   structured_grid%lye = structured_grid%lys + structured_grid%nly
   277)   structured_grid%lze = structured_grid%lzs + structured_grid%nlz
   278)   structured_grid%nlxy = structured_grid%nlx * structured_grid%nly
   279)   structured_grid%nlxz = structured_grid%nlx * structured_grid%nlz
   280)   structured_grid%nlyz = structured_grid%nly * structured_grid%nlz
   281)   structured_grid%nlmax = structured_grid%nlx * structured_grid%nly * structured_grid%nlz
   282)      
   283)   ! get ghosted corner information
   284)   call DMDAGetGhostCorners(da, structured_grid%gxs, &
   285)       structured_grid%gys, structured_grid%gzs, structured_grid%ngx, &
   286)       structured_grid%ngy, structured_grid%ngz, ierr);CHKERRQ(ierr)
   287)      
   288)   structured_grid%gxe = structured_grid%gxs + structured_grid%ngx
   289)   structured_grid%gye = structured_grid%gys + structured_grid%ngy
   290)   structured_grid%gze = structured_grid%gzs + structured_grid%ngz
   291)   structured_grid%ngxy = structured_grid%ngx * structured_grid%ngy
   292)   structured_grid%ngxz = structured_grid%ngx * structured_grid%ngz
   293)   structured_grid%ngyz = structured_grid%ngy * structured_grid%ngz
   294)   structured_grid%ngmax = structured_grid%ngx * structured_grid%ngy * structured_grid%ngz
   295)   
   296)   structured_grid%global_offset = 0
   297)   call MPI_Exscan(structured_grid%nlmax,structured_grid%global_offset, &
   298)                   ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)  
   299)    
   300) end subroutine StructGridComputeLocalBounds
   301) 
   302) ! ************************************************************************** !
   303) 
   304) subroutine StructGridCreateVecFromDM(da,vector,vector_type)
   305)   ! 
   306)   ! Creates a PETSc vector from a DM
   307)   ! 
   308)   ! Author: Glenn Hammond
   309)   ! Date: 02/08/08
   310)   ! 
   311) 
   312)   implicit none
   313) 
   314)   DM :: da
   315)   Vec :: vector
   316)   PetscInt :: vector_type
   317)   
   318)   PetscErrorCode :: ierr
   319) 
   320)   select case (vector_type)
   321)     case(GLOBAL)
   322)       call DMCreateGlobalVector(da,vector,ierr);CHKERRQ(ierr)
   323)     case(LOCAL)
   324)       call DMCreateLocalVector(da,vector,ierr);CHKERRQ(ierr)
   325)     case(NATURAL)
   326)       call DMDACreateNaturalVector(da,vector,ierr);CHKERRQ(ierr)
   327)   end select
   328) 
   329) end subroutine StructGridCreateVecFromDM
   330) 
   331) ! ************************************************************************** !
   332) 
   333) subroutine StructGridReadDXYZ(structured_grid,input,option)
   334)   ! 
   335)   ! Reads structured grid spacing from input file
   336)   ! 
   337)   ! Author: Glenn Hammond
   338)   ! Date: 10/23/07
   339)   ! 
   340) 
   341)   use Utility_module, only : UtilityReadArray
   342)   use Option_module
   343)   use Input_Aux_module
   344)   
   345)   implicit none
   346)   
   347)   type(grid_structured_type) :: structured_grid
   348)   type(input_type), pointer :: input
   349)   type(option_type) :: option
   350)   character(len=MAXSTRINGLENGTH) :: string
   351)   
   352)   PetscInt :: i
   353) 
   354)   ! the arrays are allocated within UtilityReadArray
   355)   string = 'DXYZ : X'
   356)   call UtilityReadArray(structured_grid%dx_global, &
   357)                         structured_grid%nx,string,input,option)
   358)   string = 'DXYZ : Y'
   359)   call UtilityReadArray(structured_grid%dy_global, &
   360)                         structured_grid%ny,string,input,option)
   361)   string = 'DXYZ : Z'
   362)   call UtilityReadArray(structured_grid%dz_global, &
   363)                         structured_grid%nz,string,input,option)
   364)     
   365)   if (OptionPrintToFile(option)) then
   366)     write(option%fid_out,'(/," *DXYZ ")')
   367)     write(option%fid_out,'("  dx  ",/,(1p10e12.4))') &
   368)       (structured_grid%dx_global(i),i=1,structured_grid%nx)
   369)     write(option%fid_out,'("  dy  ",/,(1p10e12.4))') &
   370)       (structured_grid%dy_global(i),i=1,structured_grid%ny)
   371)     write(option%fid_out,'("  dz  ",/,(1p10e12.4))') &
   372)       (structured_grid%dz_global(i),i=1,structured_grid%nz)
   373)   endif
   374) 
   375) end subroutine StructGridReadDXYZ
   376) 
   377) ! ************************************************************************** !
   378) 
   379) subroutine StructGridComputeSpacing(structured_grid,origin_global,option)
   380)   ! 
   381)   ! Computes structured grid spacing
   382)   ! 
   383)   ! Author: Glenn Hammond
   384)   ! Date: 10/26/07
   385)   ! 
   386) 
   387)   use Option_module
   388)   
   389)   implicit none
   390)   
   391)   type(grid_structured_type) :: structured_grid
   392)   PetscReal :: origin_global(3)
   393)   type(option_type) :: option
   394)   
   395)   PetscInt :: i, j, k, ghosted_id
   396)   PetscReal :: tempreal
   397)   PetscErrorCode :: ierr
   398) 
   399)   allocate(structured_grid%dxg_local(structured_grid%ngx))
   400)   structured_grid%dxg_local = 0.d0
   401)   allocate(structured_grid%dyg_local(structured_grid%ngy))
   402)   structured_grid%dyg_local = 0.d0
   403)   allocate(structured_grid%dzg_local(structured_grid%ngz))
   404)   structured_grid%dzg_local = 0.d0
   405)   
   406)   if (.not.associated(structured_grid%dx_global)) then
   407)     ! indicates that the grid spacings still need to be computed
   408)     if (structured_grid%bounds(1,1) < -1.d19) then 
   409)       option%io_buffer = 'Bounds have not been set for grid and DXYZ ' // & 
   410)         'does not exist'
   411)       call printErrMsg(option)
   412)     endif
   413)     allocate(structured_grid%dx_global(structured_grid%nx))
   414)     allocate(structured_grid%dy_global(structured_grid%ny))
   415)     allocate(structured_grid%dz_global(structured_grid%nz))
   416)       
   417)     select case(structured_grid%itype)
   418)       case(CARTESIAN_GRID)
   419)         structured_grid%dx_global = &
   420)           (structured_grid%bounds(X_DIRECTION,UPPER)- &
   421)            structured_grid%bounds(X_DIRECTION,LOWER)) / &
   422)           dble(structured_grid%nx)
   423)         structured_grid%dy_global = &
   424)           (structured_grid%bounds(Y_DIRECTION,UPPER)- &
   425)            structured_grid%bounds(Y_DIRECTION,LOWER)) / &
   426)           dble(structured_grid%ny)
   427)         structured_grid%dz_global = &
   428)           (structured_grid%bounds(Z_DIRECTION,UPPER)- &
   429)            structured_grid%bounds(Z_DIRECTION,LOWER)) / &
   430)           dble(structured_grid%nz)
   431)       case(CYLINDRICAL_GRID)
   432)         structured_grid%dx_global = &
   433)           (structured_grid%bounds(X_DIRECTION,UPPER)- &
   434)            structured_grid%bounds(X_DIRECTION,LOWER)) / &
   435)           dble(structured_grid%nx)
   436)         structured_grid%dy_global = 1.d0
   437)         structured_grid%dz_global = &
   438)           (structured_grid%bounds(Z_DIRECTION,UPPER)- &
   439)            structured_grid%bounds(Z_DIRECTION,LOWER)) / &
   440)           dble(structured_grid%nz)
   441)       case(SPHERICAL_GRID)
   442)         structured_grid%dx_global = &
   443)           (structured_grid%bounds(X_DIRECTION,UPPER)- &
   444)            structured_grid%bounds(X_DIRECTION,LOWER)) / &
   445)           dble(structured_grid%nx)
   446)         structured_grid%dy_global = 1.d0
   447)         structured_grid%dz_global = 1.d0
   448)     end select
   449)   else
   450)     ! x-direction
   451)     if (structured_grid%itype == CARTESIAN_GRID .or. &
   452)         structured_grid%itype == CYLINDRICAL_GRID .or. &
   453)         structured_grid%itype == SPHERICAL_GRID) then
   454)       tempreal = origin_global(X_DIRECTION)
   455)       structured_grid%bounds(X_DIRECTION,LOWER) = tempreal
   456)       do i = 1, structured_grid%nx
   457)         tempreal = tempreal + structured_grid%dx_global(i)
   458)       enddo
   459)       structured_grid%bounds(X_DIRECTION,UPPER) = tempreal
   460)     endif
   461)     ! y-direction
   462)     if (structured_grid%itype == CARTESIAN_GRID) then
   463)       tempreal = origin_global(Y_DIRECTION)
   464)       structured_grid%bounds(Y_DIRECTION,LOWER) = tempreal
   465)       do j = 1, structured_grid%ny
   466)         tempreal = tempreal + structured_grid%dy_global(j)
   467)       enddo
   468)       structured_grid%bounds(Y_DIRECTION,UPPER) = tempreal
   469)     else
   470)       structured_grid%bounds(Y_DIRECTION,LOWER) = 0.d0
   471)       structured_grid%bounds(Y_DIRECTION,UPPER) = 1.d0
   472)     endif
   473)     ! z-direction    
   474)     if (structured_grid%itype == CARTESIAN_GRID .or. &
   475)         structured_grid%itype == CYLINDRICAL_GRID) then
   476)       tempreal = origin_global(Z_DIRECTION)
   477)       structured_grid%bounds(Z_DIRECTION,LOWER) = tempreal
   478)       do k = 1, structured_grid%nz
   479)         tempreal = tempreal + structured_grid%dz_global(k)
   480)       enddo
   481)       structured_grid%bounds(Z_DIRECTION,UPPER) = tempreal
   482)     else
   483)       structured_grid%bounds(Y_DIRECTION,LOWER) = 0.d0
   484)       structured_grid%bounds(Z_DIRECTION,UPPER) = 1.d0
   485)     endif
   486)   endif
   487) 
   488)   option%io_buffer = 'Domain Bounds (x y z):'
   489)   call printMsg(option)
   490)   write(option%io_buffer,'(2x,3es18.10)') structured_grid%bounds(:,LOWER)
   491)   call printMsg(option)
   492)   write(option%io_buffer,'(2x,3es18.10)') structured_grid%bounds(:,UPPER)
   493)   call printMsg(option)
   494) 
   495)   structured_grid%dxg_local(1:structured_grid%ngx) = &
   496)     structured_grid%dx_global(structured_grid%gxs+1:structured_grid%gxe)
   497)   structured_grid%dyg_local(1:structured_grid%ngy) = &
   498)     structured_grid%dy_global(structured_grid%gys+1:structured_grid%gye)
   499)   structured_grid%dzg_local(1:structured_grid%ngz) = &
   500)     structured_grid%dz_global(structured_grid%gzs+1:structured_grid%gze)
   501)         
   502)   allocate(structured_grid%dx(structured_grid%ngmax))
   503)   structured_grid%dx = 0.d0
   504)   allocate(structured_grid%dy(structured_grid%ngmax))
   505)   structured_grid%dy = 0.d0
   506)   allocate(structured_grid%dz(structured_grid%ngmax))
   507)   structured_grid%dz = 0.d0
   508)  
   509)   do k = 1, structured_grid%ngz
   510)     do j = 1, structured_grid%ngy
   511)       do i = 1, structured_grid%ngx
   512)         ghosted_id = i+(j-1)*structured_grid%ngx+(k-1)*structured_grid%ngxy
   513)         structured_grid%dx(ghosted_id) = structured_grid%dxg_local(i)
   514)         structured_grid%dy(ghosted_id) = structured_grid%dyg_local(j)
   515)         structured_grid%dz(ghosted_id) = structured_grid%dzg_local(k)
   516)       enddo
   517)     enddo
   518)   enddo
   519)   
   520) end subroutine StructGridComputeSpacing
   521) 
   522) ! ************************************************************************** !
   523) 
   524) subroutine StructGridComputeCoord(structured_grid,option,origin_global, &
   525)                                       grid_x,grid_y,grid_z, &
   526)                                       x_min,x_max,y_min,y_max,z_min,z_max)
   527)   ! 
   528)   ! Computes structured coordinates in x,y,z
   529)   ! 
   530)   ! Author: Glenn Hammond
   531)   ! Date: 10/24/07
   532)   ! 
   533) 
   534)   use Option_module
   535)   
   536) implicit none
   537)   
   538)   type(grid_structured_type) :: structured_grid
   539)   type(option_type) :: option
   540)   PetscReal :: origin_global(3)
   541)   PetscReal :: grid_x(:), grid_y(:), grid_z(:)
   542)   PetscReal :: x_min, x_max, y_min, y_max, z_min, z_max
   543) 
   544) ! PetscErrorCode :: ierr
   545)   PetscInt :: i, j, k, ghosted_id
   546)   PetscReal :: x, y, z
   547)   PetscInt :: prevnode
   548) 
   549)   x_min = origin_global(X_DIRECTION)
   550)   y_min = origin_global(Y_DIRECTION)
   551)   z_min = origin_global(Z_DIRECTION)
   552)     
   553)   do i=1,structured_grid%lxs
   554)     x_min = x_min + structured_grid%dx_global(i)
   555)   enddo
   556)   do j=1,structured_grid%lys
   557)     y_min = y_min + structured_grid%dy_global(j)
   558)   enddo
   559)   do k=1,structured_grid%lzs
   560)     z_min = z_min + structured_grid%dz_global(k)
   561)   enddo
   562)    
   563)   ! set min and max bounds of domain in coordinate directions
   564)   structured_grid%local_origin(X_DIRECTION) = x_min
   565)   structured_grid%local_origin(Y_DIRECTION) = y_min
   566)   structured_grid%local_origin(Z_DIRECTION) = z_min
   567)   x_max = x_min
   568)   y_max = y_min
   569)   z_max = z_min
   570)   
   571)   ! there is an issue with cumulative round off error where the summed 
   572)   ! maximum extend may not match the global bound on the domain. this is 
   573)   ! an issue if a region (e.g. point is located on the upper global bound.
   574)   ! to present this, we set the maximum extend (at the global upper 
   575)   ! boundary) to the non-summed value
   576)   if (structured_grid%lxe < structured_grid%nx) then
   577)     do i=structured_grid%istart,structured_grid%iend
   578)       x_max = x_max + structured_grid%dxg_local(i+1)
   579)     enddo
   580)   else 
   581)     x_max = structured_grid%bounds(X_DIRECTION,UPPER) 
   582)   endif
   583)   if (structured_grid%lye < structured_grid%ny) then
   584)     do j=structured_grid%jstart,structured_grid%jend
   585)       y_max = y_max + structured_grid%dyg_local(j+1)
   586)     enddo
   587)   else 
   588)     y_max = structured_grid%bounds(Y_DIRECTION,UPPER) 
   589)   endif
   590)   if (structured_grid%lze < structured_grid%nz) then
   591)     do k=structured_grid%kstart,structured_grid%kend
   592)       z_max = z_max + structured_grid%dzg_local(k+1)
   593)     enddo
   594)   else 
   595)     z_max = structured_grid%bounds(Z_DIRECTION,UPPER) 
   596)   endif
   597) 
   598) ! fill in grid cell coordinates
   599)   ghosted_id = 0
   600)   if (structured_grid%kstart > 0) then
   601)     z = structured_grid%local_origin(Z_DIRECTION) - &
   602)         0.5d0*structured_grid%dzg_local(1)
   603)   else
   604)     z = structured_grid%local_origin(Z_DIRECTION) + &
   605)         0.5d0*structured_grid%dzg_local(1)
   606)   endif
   607)   do k=1, structured_grid%ngz
   608)     if (structured_grid%jstart > 0) then
   609)       y = structured_grid%local_origin(Y_DIRECTION) - &
   610)           0.5d0*structured_grid%dyg_local(1)
   611)     else
   612)       y = structured_grid%local_origin(Y_DIRECTION) + &
   613)           0.5d0*structured_grid%dyg_local(1)
   614)     endif
   615)     do j=1, structured_grid%ngy
   616)       if (structured_grid%istart > 0) then
   617)         x = structured_grid%local_origin(X_DIRECTION) - &
   618)             0.5d0*structured_grid%dxg_local(1)
   619)       else
   620)         x = structured_grid%local_origin(X_DIRECTION) + &
   621)             0.5d0*structured_grid%dxg_local(1)
   622)       endif
   623)       do i=1, structured_grid%ngx
   624)         ghosted_id = ghosted_id + 1
   625)         grid_x(ghosted_id) = x
   626)         grid_y(ghosted_id) = y
   627)         grid_z(ghosted_id) = z
   628)         if (i < structured_grid%ngx) &
   629)           x = x + &
   630)              0.5d0*(structured_grid%dxg_local(i)+structured_grid%dxg_local(i+1))
   631)       enddo
   632)       if (j < structured_grid%ngy) &
   633)         y = y + &
   634)             0.5d0*(structured_grid%dyg_local(j)+structured_grid%dyg_local(j+1))
   635)     enddo
   636)     if (k < structured_grid%ngz) &
   637)       z = z + &
   638)           0.5d0*(structured_grid%dzg_local(k)+structured_grid%dzg_local(k+1))
   639)   enddo
   640)     
   641) end subroutine StructGridComputeCoord
   642) 
   643) ! ************************************************************************** !
   644) 
   645) subroutine StructGridGetIJKFromCoordinate(structured_grid,x,y,z,i,j,k)
   646)   ! 
   647)   ! Finds local, non-ghosted i,j,k indices for
   648)   ! grid cell encompassing coordinate
   649)   ! 
   650)   ! Author: Glenn Hammond
   651)   ! Date: 08/27/08
   652)   ! 
   653) 
   654)   use Option_module
   655)   
   656)   implicit none
   657)     
   658)   type(grid_structured_type) :: structured_grid
   659)   type(option_type) :: option
   660)   PetscInt :: i, j, k
   661)   PetscInt :: i_local, j_local, k_local
   662)   PetscInt :: i_ghosted, j_ghosted, k_ghosted
   663)   PetscReal :: x, y, z
   664)   
   665)   PetscReal :: x_upper_face, y_upper_face, z_upper_face
   666) 
   667)   i = -1
   668)   j = -1
   669)   k = -1
   670) 
   671)   x_upper_face = structured_grid%local_origin(X_DIRECTION)
   672)   i_local = 1
   673)   do i_ghosted=structured_grid%istart,structured_grid%iend
   674)     if (x >= x_upper_face .and. &                   ! since i_ghosted is zero-based
   675)         x <= x_upper_face+structured_grid%dxg_local(i_ghosted+1)) then
   676)       ! test to prevent multiple procs from including a coordinate located on
   677)       ! boundary of local decomposition shared by two procs
   678)       ! if first cell in x-dir on proc
   679)       if (i_ghosted == structured_grid%istart) then
   680)         ! located on upwind boundary and ghosted
   681)           if (x == x_upper_face .and. &
   682)             structured_grid%lxs /= structured_grid%gxs) exit
   683)       endif
   684)       i = i_local 
   685)       exit
   686)     endif
   687)     i_local = i_local + 1
   688)     x_upper_face = x_upper_face + structured_grid%dxg_local(i_ghosted+1)
   689)   enddo
   690)   y_upper_face = structured_grid%local_origin(Y_DIRECTION)
   691)   j_local = 1
   692)   do j_ghosted=structured_grid%jstart,structured_grid%jend
   693)     if (y >= y_upper_face .and. &
   694)         y <= y_upper_face+structured_grid%dyg_local(j_ghosted+1)) then
   695)       ! test to prevent multiple procs from including a coordinate located on
   696)       ! boundary of local decomposition shared by two procs
   697)       ! if first cell in y-dir on proc
   698)       if (j_ghosted == structured_grid%jstart) then
   699)         ! located on upwind boundary and ghosted
   700)         if (y == y_upper_face .and. &
   701)             structured_grid%lys /= structured_grid%gys) exit
   702)       endif
   703)       j = j_local
   704)       exit
   705)     endif
   706)     j_local = j_local + 1
   707)     y_upper_face = y_upper_face + structured_grid%dyg_local(j_ghosted+1)
   708)   enddo
   709)   z_upper_face = structured_grid%local_origin(Z_DIRECTION)
   710)   k_local = 1
   711)   do k_ghosted=structured_grid%kstart,structured_grid%kend
   712)     if (z >= z_upper_face .and. &
   713)         z <= z_upper_face+structured_grid%dzg_local(k_ghosted+1)) then
   714)       ! test to prevent multiple procs from including a coordinate located on
   715)       ! boundary of local decomposition shared by two procs
   716)       ! if first cell in z-dir on proc
   717)       if (k_ghosted == structured_grid%kstart) then
   718)         ! if located on upwind boundary and ghosted, skip
   719)         if (z == z_upper_face .and. &
   720)           structured_grid%lzs /= structured_grid%gzs) exit
   721)       endif          
   722)       k = k_local
   723)       exit
   724)     endif
   725)     k_local = k_local + 1
   726)     z_upper_face = z_upper_face + structured_grid%dzg_local(k_ghosted+1)
   727)   enddo
   728)     
   729) end subroutine StructGridGetIJKFromCoordinate
   730) 
   731) ! ************************************************************************** !
   732) 
   733) subroutine StructGridGetIJKFromLocalID(structured_grid,local_id,i,j,k)
   734)   ! 
   735)   ! Finds i,j,k indices for grid cell defined by
   736)   ! local_id
   737)   ! 
   738)   ! Author: Glenn Hammond
   739)   ! Date: 04/11/08
   740)   ! 
   741) 
   742)   use Option_module
   743)   
   744)   implicit none
   745)   
   746)   type(grid_structured_type) :: structured_grid
   747)   type(option_type) :: option
   748)   PetscInt :: local_id
   749)   PetscReal :: i, j, k
   750)   
   751)   k= int((local_id-1)/structured_grid%nlxy) + 1
   752)   j= int(mod((local_id-1),structured_grid%nlxy)/structured_grid%nlx) + 1
   753)   i= mod(mod((local_id-1),structured_grid%nlxy),structured_grid%nlx) + 1  
   754)   
   755) end subroutine StructGridGetIJKFromLocalID
   756) 
   757) ! ************************************************************************** !
   758) 
   759) subroutine StructGridGetIJKFromGhostedID(structured_grid,ghosted_id,i,j,k)
   760)   ! 
   761)   ! Finds i,j,k indices for grid cell defined by
   762)   ! a ghosted id
   763)   ! 
   764)   ! Author: Glenn Hammond
   765)   ! Date: 04/11/08
   766)   ! 
   767) 
   768)   use Option_module
   769)   
   770)   implicit none
   771)   
   772)   type(grid_structured_type) :: structured_grid
   773)   type(option_type) :: option
   774)   PetscInt :: ghosted_id
   775)   PetscInt :: i, j, k
   776)   
   777)   k= int((ghosted_id-1)/structured_grid%ngxy) + 1
   778)   j= int(mod((ghosted_id-1),structured_grid%ngxy)/structured_grid%ngx) + 1
   779)   i= mod(mod((ghosted_id-1),structured_grid%ngxy),structured_grid%ngx) + 1  
   780)   
   781) end subroutine StructGridGetIJKFromGhostedID
   782) 
   783) ! ************************************************************************** !
   784) 
   785) function StructGridGetLocalIDFromIJK(structured_grid,i,j,k)
   786)   ! 
   787)   ! Finds local_id for grid cell defined by
   788)   ! i,j,k indices
   789)   ! 
   790)   ! Author: Glenn Hammond
   791)   ! Date: 01/28/11
   792)   ! 
   793) 
   794)   use Option_module
   795)   
   796)   implicit none
   797)   
   798)   type(grid_structured_type) :: structured_grid
   799)   type(option_type) :: option
   800)   PetscInt :: i, j, k
   801)   
   802)   PetscInt :: StructGridGetLocalIDFromIJK
   803)   
   804)   StructGridGetLocalIDFromIJK = &
   805)     i+(j-1)*structured_grid%nlx+(k-1)*structured_grid%nlxy
   806)   
   807) end function StructGridGetLocalIDFromIJK
   808) 
   809) ! ************************************************************************** !
   810) 
   811) function StructGridGetGhostedIDFromIJK(structured_grid,i,j,k)
   812)   ! 
   813)   ! Finds ghosted_id for grid cell defined by
   814)   ! i,j,k indices
   815)   ! 
   816)   ! Author: Glenn Hammond
   817)   ! Date: 01/28/11
   818)   ! 
   819) 
   820)   use Option_module
   821)   
   822)   implicit none
   823)   
   824)   type(grid_structured_type) :: structured_grid
   825)   type(option_type) :: option
   826)   PetscInt :: i, j, k
   827)   
   828)   PetscInt :: StructGridGetGhostedIDFromIJK
   829)   
   830)   StructGridGetGhostedIDFromIJK = &
   831)     i+(j-1)*structured_grid%ngx+(k-1)*structured_grid%ngxy
   832)   
   833) end function StructGridGetGhostedIDFromIJK
   834) 
   835) ! ************************************************************************** !
   836) 
   837) function StructGridComputeInternConnect(structured_grid, xc, yc, zc, option)
   838)   ! 
   839)   ! computes internal connectivity of a
   840)   ! structured grid
   841)   ! 
   842)   ! Author: Glenn Hammond
   843)   ! Date: 10/17/07
   844)   ! 
   845) 
   846)   use Connection_module
   847)   use Option_module
   848)   
   849)   implicit none
   850) 
   851) !  PetscReal :: radius(:)
   852)   type(connection_set_type), pointer :: StructGridComputeInternConnect
   853)   type(option_type) :: option
   854)   type(grid_structured_type) :: structured_grid
   855)   PetscReal, pointer :: xc(:),yc(:),zc(:)
   856)   
   857)   PetscReal, parameter :: Pi=3.141592653590d0
   858)   
   859)   PetscInt :: i, j, k, iconn, id_up, id_dn, id_up2, id_dn2
   860)   PetscInt :: nconn
   861)   PetscInt :: lenx, leny, lenz
   862)   PetscInt :: tvd_ghost_offset, ghost_count
   863)   PetscReal :: dist_up, dist_dn
   864)   PetscReal :: r1, r2
   865)   type(connection_set_type), pointer :: connections, connections_2
   866) 
   867)   PetscErrorCode :: ierr
   868)   PetscReal, pointer :: radius(:)
   869)   PetscInt, pointer :: int_array1(:), int_array2(:),int_array3(:),int_array4(:),int_array5(:),index(:)
   870)   PetscInt :: count
   871) 
   872)   radius => xc
   873) 
   874)   ! the adjustments in the case of AMR are based on the PIMS code adjustments by LC
   875)   nconn = (structured_grid%ngx-1)*structured_grid%nly*structured_grid%nlz+ &
   876)           structured_grid%nlx*(structured_grid%ngy-1)*structured_grid%nlz+ &
   877)           structured_grid%nlx*structured_grid%nly*(structured_grid%ngz-1)
   878)   
   879)   structured_grid%nlmax_faces = 0
   880)   structured_grid%ngmax_faces = 0
   881) 
   882)   lenx = structured_grid%ngx - 1
   883)   leny = structured_grid%ngy - 1
   884)   lenz = structured_grid%ngz - 1
   885) 
   886)   connections => ConnectionCreate(nconn,INTERNAL_CONNECTION_TYPE)
   887)   
   888)   ! if using higher order advection, allocate associated arrays
   889)   if (option%itranmode == EXPLICIT_ADVECTION .and. &
   890)       option%transport%tvd_flux_limiter /= 1) then  ! 1 = upwind
   891)     allocate(connections%id_up2(size(connections%id_up)))
   892)     allocate(connections%id_dn2(size(connections%id_dn)))
   893)     connections%id_up2 = 0
   894)     connections%id_dn2 = 0
   895)   endif
   896) 
   897)   iconn = 0
   898)   tvd_ghost_offset = 0
   899)   ghost_count = 0
   900)   ! x-connections
   901)   if (structured_grid%ngx > 1) then
   902)     select case(structured_grid%itype)
   903)       case(CARTESIAN_GRID)
   904)         do k = structured_grid%kstart, structured_grid%kend
   905)           do j = structured_grid%jstart, structured_grid%jend
   906)             do i = 1, lenx
   907)               iconn = iconn+1
   908)               id_up = i + j * structured_grid%ngx + k * structured_grid%ngxy
   909)               id_dn = id_up + 1
   910)               connections%id_up(iconn) = id_up
   911)               connections%id_dn(iconn) = id_dn
   912)               
   913)               if (associated(connections%id_up2)) then
   914)                 if (i == 1) then
   915)                   ! id_up indexes tvd_ghost_vec, see StructGridCreateTVDGhosts() 
   916) !                  id_up2 = 1 + j + k*structured_grid%nly
   917)                   ghost_count = ghost_count + 1
   918)                   id_up2 = ghost_count
   919)                   connections%id_up2(iconn) = -id_up2
   920)                 else
   921)                   connections%id_up2(iconn) = id_up - 1
   922)                 endif
   923)                 if (i == lenx) then
   924)                   ! id_dn indexes tvd_ghost_vec, see StructGridCreateTVDGhosts() 
   925) !                  id_dn2 = 1 + j + k*structured_grid%nly + structured_grid%nlyz
   926)                   id_dn2 = ghost_count + structured_grid%nlyz
   927)                   connections%id_dn2(iconn) = -id_dn2
   928)                 else
   929)                   connections%id_dn2(iconn) = id_dn + 1
   930)                 endif
   931)               endif
   932)               
   933)               connections%dist(-1:3,iconn) = 0.d0
   934)               dist_up = 0.5d0*structured_grid%dx(id_up)
   935)               dist_dn = 0.5d0*structured_grid%dx(id_dn)
   936)               connections%dist(-1,iconn) = dist_up/(dist_up+dist_dn)
   937)               connections%dist(0,iconn) = dist_up+dist_dn
   938)               connections%dist(1,iconn) = 1.d0  ! x component of unit vector
   939)               connections%area(iconn) = structured_grid%dy(id_up)* &
   940)                                         structured_grid%dz(id_up)
   941)             enddo
   942)           enddo
   943)         enddo
   944)         tvd_ghost_offset = 2*structured_grid%nlyz ! west & east
   945)         ghost_count = tvd_ghost_offset
   946)       case(CYLINDRICAL_GRID)
   947)         do k = structured_grid%kstart, structured_grid%kend
   948)           do j = structured_grid%jstart, structured_grid%jend
   949)             do i = 1, lenx
   950)               iconn = iconn+1
   951)               id_up = i + j * structured_grid%ngx + k * structured_grid%ngxy
   952)               id_dn = id_up + 1
   953)               connections%id_up(iconn) = id_up
   954)               connections%id_dn(iconn) = id_dn
   955)               connections%dist(-1:3,iconn) = 0.d0
   956)               dist_up = 0.5d0*structured_grid%dx(id_up)
   957)               dist_dn = 0.5d0*structured_grid%dx(id_dn)
   958)               connections%dist(-1,iconn) = dist_up/(dist_up+dist_dn)
   959)               connections%dist(0,iconn) = dist_up+dist_dn
   960)               connections%dist(1,iconn) = 1.d0  ! x component of unit vector
   961)               connections%area(iconn) = 2.d0 * pi * (radius(id_up)+0.5d0*structured_grid%dx(id_up))* &
   962)                                         structured_grid%dz(id_up)
   963)             enddo
   964)           enddo
   965)         enddo
   966)       case(SPHERICAL_GRID)
   967)         do k = structured_grid%kstart, structured_grid%kend
   968)           do j = structured_grid%jstart, structured_grid%jend
   969)             do i = 1, lenx
   970)               iconn = iconn+1
   971)               id_up = i + j * structured_grid%ngx + k * structured_grid%ngxy
   972)               id_dn = id_up + 1
   973)               connections%id_up(iconn) = id_up
   974)               connections%id_dn(iconn) = id_dn
   975)               connections%dist(-1:3,iconn) = 0.d0
   976)               dist_up = 0.5d0*structured_grid%dx(id_up)
   977)               dist_dn = 0.5d0*structured_grid%dx(id_dn)
   978)               connections%dist(-1,iconn) = dist_up/(dist_up+dist_dn)
   979)               connections%dist(0,iconn) = dist_up+dist_dn
   980)               connections%dist(1,iconn) = 1.d0  ! x component of unit vector
   981)               connections%area(iconn) = 4.d0 * pi * (radius(id_up)+0.5d0*structured_grid%dx(id_up))
   982)             enddo
   983)           enddo
   984)         enddo
   985)     end select
   986)   endif
   987) 
   988)   ! y-connections
   989)   if (structured_grid%ngy > 1) then
   990)     select case(structured_grid%itype)
   991)       case(CARTESIAN_GRID)
   992)         do k = structured_grid%kstart, structured_grid%kend
   993)           do i = structured_grid%istart, structured_grid%iend
   994)             do j = 1, leny
   995)               iconn = iconn+1
   996) 
   997)               id_up = i + 1 + (j-1) * structured_grid%ngx + k * structured_grid%ngxy
   998)               id_dn = id_up + structured_grid%ngx
   999)               connections%id_up(iconn) = id_up
  1000)               connections%id_dn(iconn) = id_dn
  1001)               
  1002)               if (associated(connections%id_up2)) then
  1003)                 if (j == 1) then
  1004)                   ! id_up indexes tvd_ghost_vec, see StructGridCreateTVDGhosts() 
  1005) !                  id_up2 = 1 + i + k*structured_grid%nlx + tvd_ghost_offset
  1006)                   ghost_count = ghost_count + 1
  1007)                   id_up2 = ghost_count
  1008)                   connections%id_up2(iconn) = -id_up2
  1009)                 else
  1010)                   connections%id_up2(iconn) = id_up - structured_grid%ngx
  1011)                 endif
  1012)                 if (j == leny) then
  1013)                   ! id_dn indexes tvd_ghost_vec, see StructGridCreateTVDGhosts() 
  1014) !                  id_dn2 = 1 + i + k*structured_grid%nlx + &
  1015) !                           structured_grid%nlxz + tvd_ghost_offset
  1016)                   id_dn2 = ghost_count + structured_grid%nlxz
  1017)                   connections%id_dn2(iconn) = -id_dn2
  1018)                 else
  1019)                   connections%id_dn2(iconn) = id_dn + structured_grid%ngx
  1020)                 endif
  1021)               endif
  1022)                             
  1023)               connections%dist(-1:3,iconn) = 0.d0
  1024)               dist_up = 0.5d0*structured_grid%dy(id_up)
  1025)               dist_dn = 0.5d0*structured_grid%dy(id_dn)
  1026)               connections%dist(-1,iconn) = dist_up/(dist_up+dist_dn)
  1027)               connections%dist(0,iconn) = dist_up+dist_dn
  1028)               connections%dist(2,iconn) = 1.d0  ! y component of unit vector
  1029)               connections%area(iconn) = structured_grid%dx(id_up)* &
  1030)                                     structured_grid%dz(id_up)
  1031)             enddo
  1032)           enddo
  1033)         enddo
  1034)         tvd_ghost_offset = tvd_ghost_offset + &
  1035)           2*structured_grid%nlxz ! south & north
  1036)         ghost_count = tvd_ghost_offset
  1037)       case(CYLINDRICAL_GRID)
  1038)         option%io_buffer = 'For cylindrical coordinates, NY must be equal to 1.'
  1039)         call printErrMsg(option)
  1040)       case(SPHERICAL_GRID)
  1041)         option%io_buffer = 'For spherical coordinates, NY must be equal to 1.'
  1042)         call printErrMsg(option)
  1043)     end select
  1044)   endif
  1045)       
  1046)   ! z-connections
  1047)   if (structured_grid%ngz > 1) then
  1048)     select case(structured_grid%itype)
  1049)       case(CARTESIAN_GRID)
  1050)         do j = structured_grid%jstart, structured_grid%jend
  1051)           do i = structured_grid%istart, structured_grid%iend
  1052)             do k = 1, lenz
  1053)               iconn = iconn+1
  1054) 
  1055)               id_up = i + 1 + j * structured_grid%ngx + (k-1) * &
  1056)                   structured_grid%ngxy 
  1057)               id_dn = id_up + structured_grid%ngxy
  1058)               connections%id_up(iconn) = id_up
  1059)               connections%id_dn(iconn) = id_dn
  1060)               
  1061)               if (associated(connections%id_up2)) then
  1062)                 if (k == 1) then
  1063) !                  id_up2 = 1 + i + j*structured_grid%nlx + tvd_ghost_offset
  1064)                   ghost_count = ghost_count + 1
  1065)                   id_up2 = ghost_count
  1066)                   connections%id_up2(iconn) = -id_up2
  1067)                 else
  1068)                   connections%id_up2(iconn) = id_up - structured_grid%ngxy
  1069)                 endif
  1070)                 if (k == lenz) then
  1071)                   ! id_dn indexes tvd_ghost_vec, see StructGridCreateTVDGhosts() 
  1072) !                  id_dn2 = 1 + i + j*structured_grid%nlx + &
  1073) !                           structured_grid%nlxy + tvd_ghost_offset
  1074)                   id_dn2 = ghost_count + structured_grid%nlxy
  1075)                   connections%id_dn2(iconn) = -id_dn2
  1076)                 else
  1077)                   connections%id_dn2(iconn) = id_dn + structured_grid%ngxy
  1078)                 endif
  1079)               endif
  1080)                                  
  1081)               connections%dist(-1:3,iconn) = 0.d0
  1082)               dist_up = 0.5d0*structured_grid%dz(id_up)
  1083)               dist_dn = 0.5d0*structured_grid%dz(id_dn)
  1084)               connections%dist(-1,iconn) = dist_up/(dist_up+dist_dn)
  1085)               connections%dist(0,iconn) = dist_up+dist_dn
  1086)               connections%dist(3,iconn) = 1.d0  ! z component of unit vector
  1087)               connections%area(iconn) = structured_grid%dx(id_up) * &
  1088)                                         structured_grid%dy(id_up)
  1089)             enddo
  1090)           enddo
  1091)         enddo
  1092)       case(CYLINDRICAL_GRID)
  1093)         do j = structured_grid%jstart, structured_grid%jend
  1094)           do i = structured_grid%istart, structured_grid%iend
  1095)             do k = 1, lenz
  1096)               iconn = iconn+1
  1097)               id_up = i + 1 + j * structured_grid%ngx + (k-1) * &
  1098)                   structured_grid%ngxy
  1099)               id_dn = id_up + structured_grid%ngxy
  1100)               connections%id_up(iconn) = id_up
  1101)               connections%id_dn(iconn) = id_dn
  1102)               connections%dist(-1:3,iconn) = 0.d0
  1103)               dist_up = 0.5d0*structured_grid%dz(id_up)
  1104)               dist_dn = 0.5d0*structured_grid%dz(id_dn)
  1105)               connections%dist(-1,iconn) = dist_up/(dist_up+dist_dn)
  1106)               connections%dist(0,iconn) = dist_up+dist_dn
  1107)               connections%dist(3,iconn) = 1.d0  ! z component of unit vector
  1108)               ! pi*(r2^2-r1^2)
  1109)               r2 = xc(id_up) + 0.5d0*structured_grid%dx(id_up)
  1110)               r1 = xc(id_up) - 0.5d0*structured_grid%dx(id_up)
  1111)               connections%area(iconn) = pi * dabs(r2*r2 - r1*r1)
  1112)             enddo
  1113)           enddo
  1114)         enddo
  1115)       case(SPHERICAL_GRID)
  1116)         option%io_buffer = 'For spherical coordinates, NZ must be equal to 1.'
  1117)         call printErrMsg(option)
  1118)   end select
  1119)   endif
  1120) 
  1121)   StructGridComputeInternConnect => connections
  1122) 
  1123) end function StructGridComputeInternConnect
  1124) 
  1125) ! ************************************************************************** !
  1126) 
  1127) subroutine StructGridPopulateConnection(radius,structured_grid,connection,iface, &
  1128)                                         iconn,ghosted_id,option)
  1129)   ! 
  1130)   ! Computes details of connection (area, dist, etc)
  1131)   ! 
  1132)   ! Author: Glenn Hammond
  1133)   ! Date: 11/09/07
  1134)   ! 
  1135) 
  1136)   use Connection_module
  1137)   use Option_module
  1138)   
  1139)   implicit none
  1140)  
  1141)   type(grid_structured_type) :: structured_grid
  1142)   type(connection_set_type) :: connection
  1143)   PetscInt :: iface
  1144)   PetscInt :: iconn
  1145)   PetscInt :: ghosted_id
  1146)   type(option_type) :: option
  1147)   PetscReal :: radius(:)
  1148)   
  1149)   PetscErrorCode :: ierr
  1150)   
  1151)   PetscReal, parameter :: Pi=3.141592653590d0
  1152)   
  1153)   select case(connection%itype)
  1154)     case(BOUNDARY_CONNECTION_TYPE)
  1155)     
  1156)       select case(iface)
  1157) 
  1158)         case(WEST_FACE,EAST_FACE)
  1159) 
  1160)           select case(structured_grid%itype)
  1161)             case(CARTESIAN_GRID)
  1162)               connection%dist(:,iconn) = 0.d0
  1163)               connection%dist(0,iconn) = 0.5d0*structured_grid%dx(ghosted_id)
  1164)               connection%area(iconn) = structured_grid%dy(ghosted_id)* &
  1165)                                    structured_grid%dz(ghosted_id)
  1166) 
  1167) 
  1168)               if (iface ==  WEST_FACE) then
  1169)                 connection%dist(1,iconn) = 1.d0
  1170)               else
  1171)                 connection%dist(1,iconn) = -1.d0
  1172)               endif
  1173)               if (associated(connection%id_dn2)) then
  1174)                   connection%id_dn2(iconn) = &
  1175)                     StructGetTVDGhostConnection(ghosted_id, &
  1176)                                                 structured_grid, &
  1177)                                                 iface,option)
  1178)               endif              
  1179)             case(CYLINDRICAL_GRID)
  1180)               connection%dist(:,iconn) = 0.d0
  1181)               connection%dist(0,iconn) = 0.5d0*structured_grid%dx(ghosted_id)
  1182)               if (iface ==  WEST_FACE) then
  1183)                 connection%dist(1,iconn) = 1.d0
  1184)                 connection%area(iconn) = 2.d0 * pi * (radius(ghosted_id)- &
  1185)                                       0.5d0*structured_grid%dx(ghosted_id)) * &
  1186)                                       structured_grid%dz(ghosted_id)
  1187)               else
  1188)                 connection%dist(1,iconn) = -1.d0
  1189)                 connection%area(iconn) = 2.d0 * pi * (radius(ghosted_id)+ &
  1190)                                       0.5d0*structured_grid%dx(ghosted_id)) * &
  1191)                                       structured_grid%dz(ghosted_id)
  1192)               endif
  1193)             case(SPHERICAL_GRID)
  1194)               connection%dist(:,iconn) = 0.d0
  1195)               connection%dist(0,iconn) = 0.5d0*structured_grid%dx(ghosted_id)
  1196)               if (iface ==  WEST_FACE) then
  1197)                 connection%dist(1,iconn) = 1.d0
  1198)                 connection%area(iconn) = 0.d0
  1199)               else
  1200)                 connection%dist(1,iconn) = -1.d0
  1201)                 connection%area(iconn) = 4.d0 * pi * (radius(ghosted_id)+ &
  1202)                                          0.5d0*structured_grid%dx(ghosted_id))
  1203)               endif
  1204)           end select
  1205) 
  1206)         case(SOUTH_FACE,NORTH_FACE)
  1207) 
  1208)           select case(structured_grid%itype)
  1209)             case(CARTESIAN_GRID)
  1210)               connection%dist(:,iconn) = 0.d0
  1211)               connection%dist(0,iconn) = 0.5d0*structured_grid%dy(ghosted_id)
  1212)               connection%area(iconn) = structured_grid%dx(ghosted_id)* &
  1213)                                    structured_grid%dz(ghosted_id)
  1214)               if (iface == SOUTH_FACE) then
  1215)                 connection%dist(2,iconn) = 1.d0
  1216)               else
  1217)                 connection%dist(2,iconn) = -1.d0
  1218)               endif
  1219)               if (associated(connection%id_dn2)) then
  1220)                   connection%id_dn2(iconn) = &
  1221)                     StructGetTVDGhostConnection(ghosted_id, &
  1222)                                                 structured_grid, &
  1223)                                                 iface,option)
  1224)               endif              
  1225)             case(CYLINDRICAL_GRID)
  1226)               print *, 'Cylindrical coordinates not applicable.'
  1227)               stop
  1228)             case(SPHERICAL_GRID)
  1229)               print *, 'Spherical coordinates not applicable.'
  1230)               stop
  1231)           end select
  1232) 
  1233)         case(BOTTOM_FACE,TOP_FACE)
  1234) 
  1235)           select case(structured_grid%itype)
  1236)             case(CARTESIAN_GRID)
  1237)               connection%dist(:,iconn) = 0.d0
  1238)               connection%dist(0,iconn) = 0.5d0*structured_grid%dz(ghosted_id)
  1239)               connection%area(iconn) = structured_grid%dx(ghosted_id)* &
  1240)                                    structured_grid%dy(ghosted_id)
  1241)               if (structured_grid%invert_z_axis) then
  1242)                 if (iface == TOP_FACE) then 
  1243)                   option%io_buffer = 'Need to ensure that direction of ' // &
  1244)                     'inverted z is correct in StructGridPopulateConnection()'
  1245)                   call printErrMsg(option)
  1246)                   connection%dist(3,iconn) = -1.d0
  1247)                 else
  1248)                   connection%dist(3,iconn) = 1.d0
  1249)                 endif
  1250)               else
  1251)                 if (iface == BOTTOM_FACE) then 
  1252)                   connection%dist(3,iconn) = 1.d0
  1253)                 else
  1254)                   connection%dist(3,iconn) = -1.d0
  1255)                 endif
  1256)               endif
  1257)               if (associated(connection%id_dn2)) then
  1258)                   connection%id_dn2(iconn) = &
  1259)                     StructGetTVDGhostConnection(ghosted_id, &
  1260)                                                 structured_grid, &
  1261)                                                 iface,option)
  1262)               endif              
  1263)             case(CYLINDRICAL_GRID)
  1264)               connection%dist(:,iconn) = 0.d0
  1265)               connection%dist(0,iconn) = 0.5d0*structured_grid%dz(ghosted_id)
  1266)               connection%area(iconn) = 2.d0 * pi * radius(ghosted_id) * &
  1267)                                         structured_grid%dx(ghosted_id)
  1268)               if (structured_grid%invert_z_axis) then
  1269)                 if (iface ==  TOP_FACE) then 
  1270)                   connection%dist(3,iconn) = 1.d0
  1271)                 else
  1272)                   connection%dist(3,iconn) = -1.d0
  1273)                 endif
  1274)               else
  1275)                 if (iface ==  TOP_FACE) then 
  1276)                   connection%dist(3,iconn) = -1.d0
  1277)                 else
  1278)                   connection%dist(3,iconn) = 1.d0
  1279)                 endif
  1280)               endif
  1281)             case(SPHERICAL_GRID)
  1282)               option%io_buffer = &
  1283)                 'Areas for spherical coordinates for z-axis not applicable.'
  1284)               call printErrMsg(option)
  1285)           end select
  1286)       end select
  1287)       if (connection%area(iconn) < 1.d-20) then
  1288)         write(option%io_buffer,*) connection%id_dn(iconn)
  1289)         option%io_buffer = &
  1290)           'Zero area in boundary connection at grid cell ' // &
  1291)           trim(adjustl(option%io_buffer)) // '.'
  1292)         call printErrMsg(option)
  1293)       endif
  1294)     case(INITIAL_CONNECTION_TYPE)
  1295)     case(SRC_SINK_CONNECTION_TYPE)
  1296)   end select
  1297)   
  1298) end subroutine StructGridPopulateConnection
  1299) 
  1300) ! ************************************************************************** !
  1301) 
  1302) subroutine StructGridComputeVolumes(radius,structured_grid,option,nL2G,volume)
  1303)   ! 
  1304)   ! Computes the volumes of cells in structured grid
  1305)   ! 
  1306)   ! Author: Glenn Hammond
  1307)   ! Date: 10/25/07
  1308)   ! 
  1309) 
  1310)   use Option_module
  1311)   
  1312)   implicit none
  1313) 
  1314) ! These includes are needed for VecRestoreArrayF90() - geh
  1315) #include "petsc/finclude/petscvec.h"
  1316) #include "petsc/finclude/petscvec.h90"
  1317)   
  1318)   type(grid_structured_type) :: structured_grid
  1319)   type(option_type) :: option
  1320)   PetscInt :: nL2G(:)
  1321)   PetscReal :: radius(:)
  1322)   Vec :: volume
  1323)   
  1324)   PetscReal, parameter :: Pi=3.141592653590d0
  1325)   
  1326)   PetscInt :: local_id, ghosted_id
  1327)   PetscReal, pointer :: volume_p(:)
  1328)   PetscReal :: r1, r2
  1329)   PetscErrorCode :: ierr
  1330)   
  1331)   call VecGetArrayF90(volume,volume_p, ierr);CHKERRQ(ierr)
  1332)   
  1333)   select case(structured_grid%itype)
  1334)     case(CARTESIAN_GRID)
  1335)       do local_id=1, structured_grid%nlmax
  1336)         ghosted_id = nL2G(local_id)
  1337)         volume_p(local_id) = structured_grid%dx(ghosted_id) * &
  1338)                              structured_grid%dy(ghosted_id) * &
  1339)                              structured_grid%dz(ghosted_id)
  1340)       enddo
  1341)     case(CYLINDRICAL_GRID)
  1342)       do local_id=1, structured_grid%nlmax
  1343)         ghosted_id = nL2G(local_id)
  1344)         ! volume = 2 pi r dr dz
  1345)         !       = pi * (r2-r1) * (r2+r1) * dz where dr = r2-r1 and 2 r = r2 + r1
  1346)         volume_p(local_id) = 2.d0 * pi * radius(ghosted_id) * &
  1347)                              structured_grid%dx(ghosted_id) * &
  1348)                              structured_grid%dz(ghosted_id)
  1349)       enddo
  1350)     case(SPHERICAL_GRID)
  1351)       do local_id=1, structured_grid%nlmax
  1352)         ghosted_id = nL2G(local_id)
  1353)         r2 = radius(ghosted_id) + 0.5d0*structured_grid%dx(ghosted_id)
  1354)         r1 = radius(ghosted_id) - 0.5d0*structured_grid%dx(ghosted_id)
  1355)         volume_p(local_id) = 4.d0/3.d0 * pi * structured_grid%dx(ghosted_id) &
  1356)                              * (r2*r2 + r2*r1 + r1*r1)
  1357)       enddo
  1358)   end select
  1359)   
  1360)   call VecRestoreArrayF90(volume,volume_p, ierr);CHKERRQ(ierr)
  1361)   
  1362)   if (option%print_to_screen .and. &
  1363)       option%mycommsize > 1 .and. option%mycommsize <= 16) then
  1364)     write(*,'(" rank= ",i3,", nlmax= ",i6,", nlx,y,z= ",3i4, &
  1365)       & ", nxs,e = ",2i4,", nys,e = ",2i4,", nzs,e = ",2i4)') &
  1366)       option%myrank,structured_grid%nlmax,structured_grid%nlx, &
  1367)         structured_grid%nly,structured_grid%nlz,structured_grid%lxs, &
  1368)         structured_grid%lxe,structured_grid%lys,structured_grid%lye, &
  1369)         structured_grid%lzs,structured_grid%lze
  1370) 
  1371)     write(*,'(" rank= ",i3,", ngmax= ",i6,", ngx,y,z= ",3i4, &
  1372)       & ", ngxs,e= ",2i4,", ngys,e= ",2i4,", ngzs,e= ",2i4)') &
  1373)       option%myrank,structured_grid%ngmax,structured_grid%ngx, &
  1374)         structured_grid%ngy,structured_grid%ngz,structured_grid%gxs, &
  1375)         structured_grid%gxe,structured_grid%gys,structured_grid%gye, &
  1376)         structured_grid%gzs,structured_grid%gze
  1377)   endif
  1378) 
  1379) end subroutine StructGridComputeVolumes
  1380) 
  1381) ! ************************************************************************** !
  1382) 
  1383) subroutine StructGridMapIndices(structured_grid,stencil_type, &
  1384)                                 nG2L,nL2G,nG2A,option)
  1385)   ! 
  1386)   ! maps global, local and natural indices of cells
  1387)   ! to each other
  1388)   ! 
  1389)   ! Author: Glenn Hammond
  1390)   ! Date: 10/24/07
  1391)   ! 
  1392) 
  1393)   use Option_module
  1394) 
  1395)   implicit none
  1396)   
  1397) #include "petsc/finclude/petscdm.h"
  1398) #include "petsc/finclude/petscdm.h90"
  1399) #include "petsc/finclude/petscdmda.h"  
  1400) 
  1401)   type(grid_structured_type) :: structured_grid
  1402)   PetscEnum :: stencil_type
  1403)   PetscInt, pointer :: nG2L(:), nL2G(:), nG2A(:)
  1404)   type(option_type) :: option
  1405) 
  1406)   PetscInt :: i, j, k, local_id, ghosted_id, natural_id, count1
  1407)   PetscErrorCode :: ierr
  1408)   
  1409)   allocate(nL2G(structured_grid%nlmax))
  1410)   allocate(nG2L(structured_grid%ngmax))
  1411)   allocate(nG2A(structured_grid%ngmax))
  1412) 
  1413)   structured_grid%istart = structured_grid%lxs-structured_grid%gxs
  1414)   structured_grid%jstart = structured_grid%lys-structured_grid%gys
  1415)   structured_grid%kstart = structured_grid%lzs-structured_grid%gzs
  1416)   structured_grid%iend = structured_grid%istart+structured_grid%nlx-1
  1417)   structured_grid%jend = structured_grid%jstart+structured_grid%nly-1
  1418)   structured_grid%kend = structured_grid%kstart+structured_grid%nlz-1
  1419) 
  1420)   ! Local <-> Ghosted Transformation
  1421)   nG2L = 0  ! Must initialize this to zero!
  1422)   nL2G = 0
  1423) 
  1424)   nG2A = 0
  1425) 
  1426)   local_id = 0
  1427)   do k=structured_grid%kstart,structured_grid%kend
  1428)     do j=structured_grid%jstart,structured_grid%jend
  1429)       do i=structured_grid%istart,structured_grid%iend
  1430)         local_id = local_id + 1
  1431)         ghosted_id = i+j*structured_grid%ngx+k*structured_grid%ngxy+1
  1432)         nL2G(local_id) = ghosted_id
  1433)         nG2L(ghosted_id) = local_id
  1434)       enddo
  1435)     enddo
  1436)   enddo
  1437) 
  1438)   do i=1,structured_grid%ngmax
  1439)     j = nG2L(i)
  1440)     if (j > 0) then
  1441)       k = nL2G(j)
  1442)       if (i /= k) then
  1443)         print *,'Error in ghost-local node numbering for ghost node =', i
  1444)         print *,'node_id_gtol(i) =', j
  1445)         print *,'node_id_ltog(node_id_gtol(i)) =', k
  1446)         stop
  1447)       endif
  1448)     endif
  1449)   enddo
  1450)   ! Local(non ghosted)->Natural(natural order starts from 0)
  1451) 
  1452)   ! if STAR stencil, need to set corner ghosted cells to -1
  1453)   if (stencil_type == DMDA_STENCIL_STAR) then
  1454)     !geh - set corner ghosted nodes to -1
  1455)     do k=1,structured_grid%ngz
  1456)       do j=1,structured_grid%ngy
  1457)         do i=1,structured_grid%ngx
  1458)           count1 = 0
  1459)           if (i == 1 .and. &
  1460)               abs(structured_grid%lxs-structured_grid%gxs) > 0) &
  1461)             count1 = count1 + 1
  1462)           if (i == structured_grid%ngx .and. &
  1463)               abs(structured_grid%gxe-structured_grid%lxe) > 0) &
  1464)             count1 = count1 + 1
  1465)           if (j == 1 .and. &
  1466)               abs(structured_grid%lys-structured_grid%gys) > 0) &
  1467)             count1 = count1 + 1
  1468)           if (j == structured_grid%ngy .and. &
  1469)               abs(structured_grid%gye-structured_grid%lye) > 0) &
  1470)             count1 = count1 + 1
  1471)           if (k == 1 .and. &
  1472)               abs(structured_grid%lzs-structured_grid%gzs) > 0) &
  1473)             count1 = count1 + 1
  1474)           if (k == structured_grid%ngz .and. &
  1475)               abs(structured_grid%gze-structured_grid%lze) > 0) &
  1476)             count1 = count1 + 1
  1477)           if (count1 > 1) then
  1478)             ghosted_id = i+(j-1)*structured_grid%ngx+(k-1)*structured_grid%ngxy
  1479)             nG2L(ghosted_id) = -1
  1480)           endif
  1481)         enddo
  1482)       enddo
  1483)     enddo
  1484)   endif
  1485) 
  1486)   ! local ghosted -> natural (1-based)
  1487)   ghosted_id = 0
  1488)   do k=1,structured_grid%ngz
  1489)     do j=1,structured_grid%ngy
  1490)       do i=1,structured_grid%ngx
  1491)         ghosted_id = ghosted_id + 1
  1492)         natural_id = i + structured_grid%gxs + & ! 1-based
  1493)                       (j-1+structured_grid%gys)*structured_grid%nx+ &
  1494)                       (k-1+structured_grid%gzs)*structured_grid%nxy
  1495)         nG2A(ghosted_id) = natural_id
  1496)       enddo
  1497)     enddo
  1498)   enddo
  1499) 
  1500) end subroutine StructGridMapIndices
  1501) 
  1502) ! ************************************************************************** !
  1503) 
  1504) subroutine StructGridGetGhostedNeighbors(structured_grid,ghosted_id, &
  1505)                                          stencil_type, &
  1506)                                          stencil_width_i,stencil_width_j, &
  1507)                                          stencil_width_k,x_count,y_count, &
  1508)                                          z_count,ghosted_neighbors, &
  1509)                                          option)
  1510)   ! 
  1511)   ! Returns an array of neighboring cells
  1512)   ! 
  1513)   ! Author: Glenn Hammond
  1514)   ! Date: 01/28/11
  1515)   ! 
  1516) 
  1517)   use Option_module
  1518) 
  1519)   implicit none
  1520) #include "petsc/finclude/petscdmda.h"
  1521)   
  1522)   type(grid_structured_type) :: structured_grid
  1523)   type(option_type) :: option
  1524)   PetscInt :: ghosted_id
  1525)   PetscEnum :: stencil_type
  1526)   PetscInt :: stencil_width_i
  1527)   PetscInt :: stencil_width_j
  1528)   PetscInt :: stencil_width_k
  1529)   PetscInt :: x_count
  1530)   PetscInt :: y_count
  1531)   PetscInt :: z_count
  1532)   PetscInt :: ghosted_neighbors(*)
  1533) 
  1534)   PetscInt :: i, j, k
  1535)   PetscInt :: icount
  1536)   PetscInt :: ii, jj, kk
  1537) 
  1538)   call StructGridGetIJKFromGhostedID(structured_grid,ghosted_id,i,j,k)
  1539)   
  1540)   x_count = 0
  1541)   y_count = 0
  1542)   z_count = 0
  1543)   icount = 0
  1544)   select case(stencil_type)
  1545)     case(DMDA_STENCIL_STAR)
  1546)       do ii = max(i-stencil_width_i,1), min(i+stencil_width_i,structured_grid%ngx)
  1547)         if (ii /= i) then
  1548)           icount = icount + 1
  1549)           x_count = x_count + 1
  1550)           ghosted_neighbors(icount) = &
  1551)             StructGridGetGhostedIDFromIJK(structured_grid,ii,j,k)
  1552)         endif
  1553)       enddo
  1554)       do jj = max(j-stencil_width_j,1), min(j+stencil_width_j,structured_grid%ngy)
  1555)         if (jj /= j) then
  1556)           icount = icount + 1
  1557)           y_count = y_count + 1
  1558)           ghosted_neighbors(icount) = &
  1559)             StructGridGetGhostedIDFromIJK(structured_grid,i,jj,k)
  1560)         endif
  1561)       enddo
  1562)       do kk = max(k-stencil_width_k,1), min(k+stencil_width_k,structured_grid%ngz)
  1563)         if (kk /= k) then
  1564)           icount = icount + 1
  1565)           z_count = z_count + 1
  1566)           ghosted_neighbors(icount) = &
  1567)             StructGridGetGhostedIDFromIJK(structured_grid,i,j,kk)
  1568)         endif
  1569)       enddo
  1570)     case(DMDA_STENCIL_BOX)
  1571)       option%io_buffer = 'DMDA_STENCIL_BOX not yet supported in ' // &
  1572)         'StructGridGetNeighbors.'
  1573)       call printErrMsg(option)
  1574)   end select
  1575) 
  1576) end subroutine StructGridGetGhostedNeighbors
  1577) 
  1578) ! ************************************************************************** !
  1579) 
  1580) subroutine StructGridGetGhostedNeighborsCorners(structured_grid,ghosted_id, &
  1581)                                          stencil_type, &
  1582)                                          stencil_width_i,stencil_width_j, &
  1583)                                          stencil_width_k, icount, &
  1584)                                          ghosted_neighbors, &
  1585)                                          option)
  1586)   ! 
  1587)   ! Returns an array of neighboring cells
  1588)   ! including the corner nodes
  1589)   ! Note that the previous subroutine does not return the corner nodes
  1590)   ! 
  1591)   ! Author: Satish Karra, LANL
  1592)   ! Date: 02/19/12
  1593)   ! 
  1594) 
  1595)   use Option_module
  1596) 
  1597)   implicit none
  1598) #include "petsc/finclude/petscdmda.h"
  1599)   
  1600)   type(grid_structured_type) :: structured_grid
  1601)   type(option_type) :: option
  1602)   PetscInt :: ghosted_id
  1603)   PetscEnum :: stencil_type
  1604)   PetscInt :: stencil_width_i
  1605)   PetscInt :: stencil_width_j
  1606)   PetscInt :: stencil_width_k
  1607)   PetscInt :: ghosted_neighbors(*)
  1608) 
  1609)   PetscInt :: i, j, k
  1610)   PetscInt :: icount
  1611)   PetscInt :: ii, jj, kk
  1612) 
  1613)   call StructGridGetIJKFromGhostedID(structured_grid,ghosted_id,i,j,k)
  1614) 
  1615)   icount = 0
  1616)   
  1617)   ! gb:08/08/13 Dependence on stencil_type is not necessary.
  1618)   !select case(stencil_type)
  1619)   !  case(DMDA_STENCIL_STAR)
  1620)       do kk = max(k-stencil_width_k,1), &
  1621)                 min(k+stencil_width_k,structured_grid%ngz)
  1622)         do jj = max(j-stencil_width_j,1), &
  1623)                   min(j+stencil_width_j,structured_grid%ngy)
  1624)           do ii = max(i-stencil_width_i,1), &
  1625)                     min(i+stencil_width_i,structured_grid%ngx)
  1626)             if (ii == i .and. jj == j .and. kk == k) then
  1627)             ! do nothing
  1628)             else
  1629)               icount = icount + 1
  1630)               ghosted_neighbors(icount) = &
  1631)               StructGridGetGhostedIDFromIJK(structured_grid,ii,jj,kk)
  1632)             endif
  1633)           enddo
  1634)         enddo          
  1635)       enddo
  1636)   !  case(DMDA_STENCIL_BOX)
  1637)   !    option%io_buffer = 'DMDA_STENCIL_BOX not yet supported in ' // &
  1638)   !      'StructGridGetNeighbors.'
  1639)   !    call printErrMsg(option)
  1640)   !end select
  1641) 
  1642) end subroutine StructGridGetGhostedNeighborsCorners
  1643) 
  1644) ! ************************************************************************** !
  1645) 
  1646) subroutine StructGridDestroy(structured_grid)
  1647)   ! 
  1648)   ! Deallocates a structured grid
  1649)   ! 
  1650)   ! Author: Glenn Hammond
  1651)   ! Date: 11/01/07
  1652)   ! 
  1653)   use Utility_module, only : DeallocateArray
  1654)   
  1655)   implicit none
  1656)   
  1657)   type(grid_structured_type), pointer :: structured_grid
  1658)   
  1659)   PetscErrorCode :: ierr
  1660)     
  1661)   if (.not.associated(structured_grid)) return
  1662)   
  1663)   call DeallocateArray(structured_grid%dx_global)
  1664)   call DeallocateArray(structured_grid%dy_global)
  1665)   call DeallocateArray(structured_grid%dz_global)
  1666)   
  1667)   call DeallocateArray(structured_grid%dxg_local)
  1668)   call DeallocateArray(structured_grid%dyg_local)
  1669)   call DeallocateArray(structured_grid%dzg_local)
  1670)   
  1671)   call DeallocateArray(structured_grid%dx)
  1672)   call DeallocateArray(structured_grid%dy)
  1673)   call DeallocateArray(structured_grid%dz)
  1674)   
  1675)   call DeallocateArray(structured_grid%cell_neighbors)
  1676)   
  1677)   deallocate(structured_grid)
  1678)   nullify(structured_grid)
  1679) 
  1680) end subroutine StructGridDestroy
  1681) 
  1682) ! ************************************************************************** !
  1683) 
  1684) subroutine StructGridCreateTVDGhosts(structured_grid,ndof,global_vec, &
  1685)                                      dm_1dof, &
  1686)                                      ghost_vec,scatter_ctx,option)
  1687)   ! 
  1688)   ! Calculates the TVD ghost vector and the
  1689)   ! associated scatter context
  1690)   ! 
  1691)   ! Author: Glenn Hammond
  1692)   ! Date: 01/28/11
  1693)   ! 
  1694) 
  1695)   use Option_module
  1696) 
  1697)   implicit none
  1698)   
  1699) #include "petsc/finclude/petscvec.h"
  1700) #include "petsc/finclude/petscvec.h90"
  1701) #include "petsc/finclude/petscviewer.h"
  1702) 
  1703)   type(grid_structured_type) :: structured_grid
  1704)   PetscInt :: ndof
  1705)   DM :: dm_1dof
  1706)   Vec :: global_vec
  1707)   Vec :: ghost_vec
  1708)   VecScatter :: scatter_ctx
  1709)   type(option_type) :: option
  1710) 
  1711)   PetscInt :: vector_size
  1712)   IS :: is_petsc
  1713)   IS :: is_ghost
  1714)   PetscInt :: icount, index, offset
  1715)   PetscInt :: increment
  1716)   PetscInt :: i, j, k
  1717)   PetscInt, allocatable :: global_indices_of_local_ghosted(:)
  1718)   PetscInt, allocatable :: global_indices_from(:)
  1719)   PetscInt, allocatable :: tvd_ghost_indices_to(:)
  1720)   ISLocalToGlobalMapping :: mapping_ltog
  1721)   PetscViewer :: viewer
  1722)   PetscErrorCode :: ierr
  1723)   
  1724)   ! structured grid has 6 sides to it
  1725)   vector_size = 0
  1726)   ! east-west
  1727)   if (structured_grid%nx > 1) then
  1728)     increment = structured_grid%nly*structured_grid%nlz
  1729)     vector_size = vector_size + 2*increment
  1730)   endif
  1731)   ! north-south
  1732)   if (structured_grid%ny > 1) then
  1733)     increment = structured_grid%nlx*structured_grid%nlz
  1734)     vector_size = vector_size + 2*increment
  1735)   endif
  1736)   ! top-bottom
  1737)   if (structured_grid%nz > 1) then
  1738)     increment = structured_grid%nlx*structured_grid%nly
  1739)     vector_size = vector_size + 2*increment
  1740)   endif
  1741)   
  1742)   if (vector_size == 0) then
  1743)     option%io_buffer = 'TVD does not handle a single grid cell.'
  1744)     call printErrMsg(option)
  1745)   endif
  1746)   
  1747)   call VecCreateSeq(PETSC_COMM_SELF,vector_size*ndof,ghost_vec, &
  1748)                     ierr);CHKERRQ(ierr)
  1749)   call VecSetBlockSize(ghost_vec,ndof,ierr);CHKERRQ(ierr)
  1750)   
  1751)   ! Create an IS composed of the petsc indexing of the ghost cells
  1752)   allocate(global_indices_from(vector_size))
  1753)   global_indices_from = UNINITIALIZED_INTEGER ! to catch bugs
  1754)   allocate(tvd_ghost_indices_to(vector_size))
  1755)   do i = 1, vector_size
  1756)     tvd_ghost_indices_to(i) = i-1
  1757)   enddo
  1758)   
  1759)   ! GEH: I'm going to play a trick here.  If I know the global index
  1760)   ! of my ghost cells, I can calculate the global index of the next
  1761)   ! layer of ghost cells in each direction since the block are 
  1762)   ! consistent through each dimension of the grid
  1763)   allocate(global_indices_of_local_ghosted(structured_grid%ngmax))
  1764)   do i = 1, structured_grid%ngmax
  1765)     global_indices_of_local_ghosted(i) = i-1
  1766)   enddo
  1767)   call DMGetLocalToGlobalMapping(dm_1dof,mapping_ltog,ierr);CHKERRQ(ierr)
  1768)   ! in and out integer arrays can be the same
  1769)   call ISLocalToGlobalMappingApply(mapping_ltog,structured_grid%ngmax, &
  1770)                                    global_indices_of_local_ghosted, &
  1771)                                    global_indices_of_local_ghosted, &
  1772)                                    ierr);CHKERRQ(ierr)
  1773)   ! leave global_indices_of_local_ghosted() in zero-based for the below
  1774)   
  1775)   ! Need to make a list of all indices that will receive updates through
  1776)   ! scatter/gather operation. Ghost cells representing physical boundaries
  1777)   ! do not need such an update.
  1778)   icount = 0
  1779) 
  1780)   if (structured_grid%nx > 1) then
  1781)     ! west
  1782)     offset = 0
  1783)     if (structured_grid%lxs /= structured_grid%gxs) offset = -1
  1784)     i = structured_grid%istart
  1785)     do k = structured_grid%kstart, structured_grid%kend
  1786)       do j = structured_grid%jstart, structured_grid%jend
  1787)         icount = icount + 1
  1788)         index = i + j*structured_grid%ngx + k*structured_grid%ngxy + 1
  1789)         global_indices_from(icount) = &
  1790)           global_indices_of_local_ghosted(index) + offset
  1791)       enddo
  1792)     enddo
  1793) 
  1794)     ! east
  1795)     offset = 0
  1796)     if (structured_grid%lxe /= structured_grid%gxe) offset = 1
  1797)     i = structured_grid%iend
  1798)     do k = structured_grid%kstart, structured_grid%kend
  1799)       do j = structured_grid%jstart, structured_grid%jend
  1800)         icount = icount + 1
  1801)         index = i + j*structured_grid%ngx + k*structured_grid%ngxy + 1
  1802)         global_indices_from(icount) = &
  1803)           global_indices_of_local_ghosted(index) + offset
  1804)       enddo
  1805)     enddo
  1806)   endif
  1807) 
  1808)   if (structured_grid%ny > 1) then
  1809)     ! south
  1810)     offset = 0
  1811)     if (structured_grid%lys /= structured_grid%gys) offset = -structured_grid%ngx
  1812)     j = structured_grid%jstart
  1813)     do k = structured_grid%kstart, structured_grid%kend
  1814)       do i = structured_grid%istart, structured_grid%iend
  1815)         icount = icount + 1
  1816)         index = i + j*structured_grid%ngx + k*structured_grid%ngxy + 1
  1817)         global_indices_from(icount) = &
  1818)           global_indices_of_local_ghosted(index) + offset
  1819)       enddo
  1820)     enddo
  1821)   
  1822)     ! north
  1823)     offset = 0
  1824)     if (structured_grid%lye /= structured_grid%gye) offset = structured_grid%ngx
  1825)     j = structured_grid%jend
  1826)     do k = structured_grid%kstart, structured_grid%kend
  1827)       do i = structured_grid%istart, structured_grid%iend
  1828)         icount = icount + 1
  1829)         index = i + j*structured_grid%ngx + k*structured_grid%ngxy + 1
  1830)         global_indices_from(icount) = &
  1831)           global_indices_of_local_ghosted(index) + offset
  1832)       enddo
  1833)     enddo
  1834)   endif
  1835)   
  1836)   if (structured_grid%nz > 1) then
  1837)     ! bottom
  1838)     offset = 0
  1839)     if (structured_grid%lzs /= structured_grid%gzs) offset = -structured_grid%ngxy
  1840)     k = structured_grid%kstart
  1841)     do j = structured_grid%jstart, structured_grid%jend
  1842)       do i = structured_grid%istart, structured_grid%iend
  1843)         icount = icount + 1
  1844)         index = i + j*structured_grid%ngx + k*structured_grid%ngxy + 1
  1845)         global_indices_from(icount) = &
  1846)           global_indices_of_local_ghosted(index) + offset
  1847)       enddo
  1848)     enddo
  1849)   
  1850)     ! top
  1851)     offset = 0
  1852)     if (structured_grid%lze /= structured_grid%gze) offset = structured_grid%ngxy
  1853)     k = structured_grid%kend
  1854)     do j = structured_grid%jstart, structured_grid%jend
  1855)       do i = structured_grid%istart, structured_grid%iend
  1856)         icount = icount + 1
  1857)         index = i + j*structured_grid%ngx + k*structured_grid%ngxy + 1
  1858)         global_indices_from(icount) = &
  1859)           global_indices_of_local_ghosted(index) + offset
  1860)       enddo
  1861)     enddo
  1862)   endif
  1863)   
  1864)   deallocate(global_indices_of_local_ghosted)
  1865) 
  1866)   if (vector_size /= icount) then
  1867)     option%io_buffer = 'Mis-count in TVD ghosting.'
  1868)     call printErrMsgByRank(option)
  1869)   endif
  1870) 
  1871)   ! since global_indices_from was base-zero, global_indices_from is base-zero.
  1872)   call ISCreateBlock(option%mycomm,ndof,vector_size, &
  1873)                       global_indices_from,PETSC_COPY_VALUES,is_petsc, &
  1874)                      ierr);CHKERRQ(ierr)
  1875)   deallocate(global_indices_from)
  1876) 
  1877) #if TVD_DEBUG
  1878)   call PetscViewerASCIIOpen(option%mycomm,'is_petsc_tvd.out', &
  1879)                             viewer,ierr);CHKERRQ(ierr)
  1880)   call ISView(is_petsc,viewer,ierr);CHKERRQ(ierr)
  1881)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
  1882) #endif
  1883) 
  1884)   ! already zero-based
  1885)   call ISCreateBlock(option%mycomm,ndof,vector_size, &
  1886)                       tvd_ghost_indices_to,PETSC_COPY_VALUES,is_ghost, &
  1887)                      ierr);CHKERRQ(ierr)
  1888)   deallocate(tvd_ghost_indices_to)
  1889) 
  1890) #if TVD_DEBUG
  1891)   call PetscViewerASCIIOpen(option%mycomm,'is_ghost_tvd.out', &
  1892)                             viewer,ierr);CHKERRQ(ierr)
  1893)   call ISView(is_ghost,viewer,ierr);CHKERRQ(ierr)
  1894)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
  1895) #endif
  1896) 
  1897)   call VecScatterCreate(global_vec,is_petsc,ghost_vec,is_ghost, &
  1898)                         scatter_ctx,ierr);CHKERRQ(ierr)
  1899) 
  1900) #if TVD_DEBUG
  1901)   call PetscViewerASCIIOpen(option%mycomm,'tvd_ghost_scatter.out',viewer, &
  1902)                             ierr);CHKERRQ(ierr)
  1903)   call VecScatterView(scatter_ctx,viewer,ierr);CHKERRQ(ierr)
  1904)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
  1905) #endif
  1906) 
  1907) end subroutine StructGridCreateTVDGhosts  
  1908) 
  1909) ! ************************************************************************** !
  1910) 
  1911) function StructGetTVDGhostConnection(ghosted_id,structured_grid,iface,option)
  1912)   ! 
  1913)   ! Returns id of tvd ghost cell for connection
  1914)   ! 
  1915)   ! Author: Glenn Hammond
  1916)   ! Date: 02/10/11
  1917)   ! 
  1918) 
  1919)   use Option_module
  1920) 
  1921)   implicit none
  1922)   
  1923)   PetscInt :: ghosted_id
  1924)   type(grid_structured_type) :: structured_grid
  1925)   PetscInt :: iface
  1926)   type(option_type) :: option
  1927)   
  1928)   PetscInt :: StructGetTVDGhostConnection
  1929)   
  1930)   character(len=MAXSTRINGLENGTH) :: string
  1931)   PetscInt :: index
  1932)   PetscInt :: offset
  1933)   PetscInt :: i, j, k
  1934)   PetscBool :: error
  1935)   
  1936)   error = PETSC_FALSE
  1937) 
  1938)   select case(iface)
  1939)     case(WEST_FACE,EAST_FACE)
  1940)       if (structured_grid%ngx > 1) then
  1941)         if (iface == WEST_FACE) then
  1942)           StructGetTVDGhostConnection = ghosted_id + 1
  1943)         else
  1944)           StructGetTVDGhostConnection = ghosted_id - 1
  1945)         endif
  1946)         return
  1947)       elseif (structured_grid%nx == 1) then
  1948)         option%io_buffer = 'Boundary condition cannot be assigned in X ' // &
  1949)           'dimension with nx = 1 with TVD.'
  1950)         error = PETSC_TRUE
  1951)       endif
  1952)     case(SOUTH_FACE,NORTH_FACE)
  1953)       if (structured_grid%ngy > 1) then
  1954)         if (iface == SOUTH_FACE) then
  1955)           StructGetTVDGhostConnection = ghosted_id + structured_grid%ngx
  1956)         else
  1957)           StructGetTVDGhostConnection = ghosted_id - structured_grid%ngx
  1958)         endif
  1959)         return
  1960)       elseif (structured_grid%ny == 1) then
  1961)         option%io_buffer = 'Boundary condition cannot be assigned in Y ' // &
  1962)           'dimension with ny = 1 with TVD.'
  1963)         error = PETSC_TRUE
  1964)       endif
  1965)     case(BOTTOM_FACE,TOP_FACE)
  1966)       if (structured_grid%ngz > 1) then
  1967)         if (iface == BOTTOM_FACE) then
  1968)           StructGetTVDGhostConnection = ghosted_id + structured_grid%ngxy
  1969)         else
  1970)           StructGetTVDGhostConnection = ghosted_id - structured_grid%ngxy
  1971)         endif
  1972)         return
  1973)       elseif (structured_grid%nz == 1) then
  1974)         option%io_buffer = 'Boundary condition cannot be assigned in Z ' // &
  1975)           'dimension with nz = 1 with TVD.'
  1976)         error = PETSC_TRUE
  1977)       endif
  1978)   end select
  1979) 
  1980)   if (error) call printErrMsg(option)
  1981)   
  1982)   call StructGridGetIJKFromGhostedID(structured_grid,ghosted_id,i,j,k)
  1983)   offset = 0
  1984)   select case(iface)
  1985)     case(WEST_FACE)
  1986)       ! ensure that connection is on boundary face
  1987)       if (i /= 1) then
  1988)         error = PETSC_TRUE
  1989)         string = 'WEST'
  1990)       endif                       ! this must be in local dimension nly
  1991)       index = j + k*structured_grid%nly
  1992)     case(EAST_FACE)
  1993)       if (i /= structured_grid%ngx) then
  1994)         error = PETSC_TRUE
  1995)         string = 'EAST'
  1996)       endif
  1997)       index = j + k*structured_grid%nly + structured_grid%nlyz
  1998)     case(SOUTH_FACE)
  1999)       if (j /= 1) then
  2000)         error = PETSC_TRUE
  2001)         string = 'SOUTH'
  2002)       endif
  2003)       if (structured_grid%nx > 1) then
  2004)         offset = 2*structured_grid%nlyz
  2005)       endif
  2006)       index = i + k*structured_grid%nlx + offset
  2007)     case(NORTH_FACE)
  2008)       if (j /= structured_grid%ngy) then
  2009)         error = PETSC_TRUE
  2010)         string = 'NORTH'
  2011)       endif
  2012)       if (structured_grid%nx > 1) then
  2013)         offset = 2*structured_grid%nlyz
  2014)       endif
  2015)       index = i + k*structured_grid%nlx + structured_grid%nlxz + offset
  2016)     case(BOTTOM_FACE)
  2017)       if (k /= 1) then
  2018)         error = PETSC_TRUE
  2019)         string = 'BOTTOM'
  2020)       endif
  2021)       if (structured_grid%nx > 1) then
  2022)         offset = 2*structured_grid%nlyz
  2023)       endif
  2024)       if (structured_grid%ny > 1) then
  2025)         offset = offset + 2*structured_grid%nlxz
  2026)       endif
  2027)       index = i + j*structured_grid%nlx + offset
  2028)     case(TOP_FACE)
  2029)       if (k /= structured_grid%ngz) then
  2030)         error = PETSC_TRUE
  2031)         string = 'TOP'
  2032)       endif
  2033)       if (structured_grid%nx > 1) then
  2034)         offset = 2*structured_grid%nlyz
  2035)       endif
  2036)       if (structured_grid%ny > 1) then
  2037)         offset = offset + 2*structured_grid%nlxz
  2038)       endif
  2039)       index = i + j*structured_grid%nlx + structured_grid%nlxy + offset
  2040)   end select
  2041)   
  2042)   if (error) then
  2043)     write(option%io_buffer, '(''StructGetTVDGhostConnection not on '', a, &
  2044)     & ''face for cell:'',3i6)') trim(string), i,j,k
  2045)     call printErrMsgByRank(option)
  2046)   endif
  2047)   
  2048)   StructGetTVDGhostConnection = -index
  2049) 
  2050) end function StructGetTVDGhostConnection
  2051) 
  2052) end module Grid_Structured_module

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody