discretization.F90       coverage:  56.67 %func     44.53 %block


     1) module Discretization_module
     2) 
     3)   use Grid_module
     4)   use Grid_Structured_module
     5)   use Grid_Unstructured_module
     6)   use Grid_Unstructured_Aux_module
     7)   use Grid_Unstructured_Explicit_module
     8)   use Grid_Unstructured_Polyhedra_module
     9)   use DM_Kludge_module
    10) 
    11)   use PFLOTRAN_Constants_module
    12) 
    13)   implicit none
    14) 
    15)   private
    16)  
    17) #include "petsc/finclude/petscsys.h"
    18) 
    19) #include "petsc/finclude/petscvec.h"
    20) #include "petsc/finclude/petscvec.h90"
    21) #include "petsc/finclude/petscmat.h"
    22) #include "petsc/finclude/petscmat.h90"
    23) #include "petsc/finclude/petscdm.h"
    24) #include "petsc/finclude/petscdm.h90"
    25) #include "petsc/finclude/petscdmda.h"
    26) #include "petsc/finclude/petscdmshell.h90"
    27) 
    28)   type, public :: discretization_type
    29)     PetscInt :: itype  ! type of discretization (e.g. structured, unstructured, etc.)
    30)     !geh: note that differentiating between implicit and explicit unstructured 
    31)     !     grids is handled within the grid%itype variable, not discritization%itype
    32)     character(len=MAXWORDLENGTH) :: ctype
    33)     PetscReal :: origin_global(3) ! origin of global domain
    34)     type(grid_type), pointer :: grid  ! pointer to a grid object
    35)     character(len=MAXSTRINGLENGTH) :: filename
    36) 
    37)     type(dm_ptr_type), pointer :: dmc_nflowdof(:), dmc_ntrandof(:)
    38)       ! Arrays containing hierarchy of coarsened DMs, for use with Galerkin 
    39)       ! multigrid.  Element i of each array is a *finer* DM than element i-1.
    40)     PetscInt :: dm_index_to_ndof(5) ! mapping between a dm_ptr to the number of degrees of freedom
    41)     type(dm_ptr_type), pointer :: dm_1dof
    42)     type(dm_ptr_type), pointer :: dm_nflowdof
    43)     type(dm_ptr_type), pointer :: dm_ntrandof
    44)     type(dm_ptr_type), pointer :: dm_n_stress_strain_dof
    45)     VecScatter :: tvd_ghost_scatter
    46)     
    47)     PetscInt :: stencil_width
    48)     PetscEnum :: stencil_type
    49)     
    50)   end type discretization_type
    51) 
    52)   public :: DiscretizationCreate, &
    53)             DiscretizationDestroy, &
    54)             DiscretizationReadRequiredCards, &
    55)             DiscretizationRead, &
    56)             DiscretizationCreateVector, &
    57)             DiscretizationDuplicateVector, &         
    58)             DiscretizationCreateJacobian, &
    59)             DiscretizationCreateInterpolation, &
    60)             DiscretizationCreateColoring, &
    61)             DiscretizationGlobalToLocal, &
    62)             DiscretizationLocalToGlobal, &
    63)             DiscretizationLocalToLocal, &
    64)             DiscretizationGlobalToNatural, &
    65)             DiscretizationNaturalToGlobal, &
    66)             DiscretizationGlobalToLocalBegin, &
    67)             DiscretizationGlobalToLocalEnd, &
    68)             DiscretizationLocalToLocalBegin, &
    69)             DiscretizationLocalToLocalEnd, &
    70)             DiscretizGlobalToNaturalBegin, &
    71)             DiscretizGlobalToNaturalEnd, &
    72)             DiscretizNaturalToGlobalBegin, &
    73)             DiscretizNaturalToGlobalEnd, &
    74)             DiscretizationCreateDMs,&
    75)             DiscretizationGetDMPtrFromIndex, &
    76)             DiscretizationUpdateTVDGhosts, &
    77)             DiscretAOApplicationToPetsc, &
    78)             DiscretizationInputRecord, &
    79)             DiscretizationPrintInfo
    80)   
    81) contains
    82) 
    83) ! ************************************************************************** !
    84) 
    85) function DiscretizationCreate()
    86)   ! 
    87)   ! Creates a structured or unstructured discretization
    88)   ! 
    89)   ! Author: Glenn Hammond
    90)   ! Date: 10/23/07
    91)   ! 
    92) 
    93)   implicit none
    94)   
    95)   type(discretization_type), pointer :: DiscretizationCreate
    96)   
    97)   type(discretization_type), pointer :: discretization
    98)   
    99)   allocate(discretization)
   100)   discretization%ctype = ''
   101)   discretization%itype = 0
   102)   discretization%origin_global = 0.d0
   103)   discretization%filename = ''
   104) 
   105)   ! nullify DM pointers
   106)   nullify(discretization%dmc_nflowdof)
   107)   nullify(discretization%dmc_ntrandof)
   108)   allocate(discretization%dm_1dof)
   109)   allocate(discretization%dm_nflowdof)
   110)   allocate(discretization%dm_ntrandof)
   111)   allocate(discretization%dm_n_stress_strain_dof)
   112)   discretization%dm_1dof%dm = 0
   113)   discretization%dm_nflowdof%dm = 0
   114)   discretization%dm_ntrandof%dm = 0
   115)   discretization%dm_n_stress_strain_dof%dm = 0
   116)   nullify(discretization%dm_1dof%ugdm)
   117)   nullify(discretization%dm_nflowdof%ugdm)
   118)   nullify(discretization%dm_ntrandof%ugdm)
   119)   nullify(discretization%dm_n_stress_strain_dof%ugdm)
   120)   
   121)   nullify(discretization%grid)
   122)   
   123)   discretization%stencil_width = 1
   124)   discretization%stencil_type = DMDA_STENCIL_STAR
   125) 
   126)   discretization%tvd_ghost_scatter = 0
   127)   
   128)   DiscretizationCreate => discretization
   129) 
   130) end function DiscretizationCreate
   131) 
   132) ! ************************************************************************** !
   133) 
   134) subroutine DiscretizationReadRequiredCards(discretization,input,option)
   135)   ! 
   136)   ! Reads a discretization from the input file
   137)   ! 
   138)   ! Author: Glenn Hammond
   139)   ! Date: 11/01/07
   140)   ! 
   141) 
   142)   use Option_module
   143)   use Input_Aux_module
   144)   use String_module
   145) 
   146)   implicit none
   147) 
   148)   type(option_type), pointer :: option
   149)   type(input_type), pointer :: input
   150)   type(discretization_type),pointer :: discretization
   151)   character(len=MAXWORDLENGTH) :: word
   152)   type(grid_type), pointer :: grid, grid2
   153)   type(grid_structured_type), pointer :: str_grid
   154)   type(grid_unstructured_type), pointer :: un_str_grid
   155)   character(len=MAXWORDLENGTH) :: structured_grid_ctype
   156)   character(len=MAXWORDLENGTH) :: unstructured_grid_ctype
   157) 
   158)   character(len=MAXSTRINGLENGTH) :: string
   159) 
   160)   PetscInt :: structured_grid_itype
   161)   PetscInt :: unstructured_grid_itype
   162)   PetscInt :: nx, ny, nz
   163)   PetscInt :: i
   164)   PetscReal :: tempreal
   165) 
   166)   nx = 0
   167)   ny = 0
   168)   nz = 0
   169) 
   170) ! we initialize the word to blanks to avoid error reported by valgrind
   171)   word = ''
   172) 
   173)   do
   174)   
   175)     call InputReadPflotranString(input,option)
   176)     if (input%ierr /= 0) exit
   177) 
   178)     if (InputCheckExit(input,option)) exit
   179) 
   180)     call InputReadWord(input,option,word,PETSC_TRUE)
   181)     call InputErrorMsg(input,option,'keyword','GRID')
   182)     call StringToUpper(word)
   183)       
   184)     select case(trim(word))
   185)       case('TYPE')
   186)         call InputReadWord(input,option,discretization%ctype,PETSC_TRUE)
   187)         call InputErrorMsg(input,option,'type','GRID')   
   188)         call StringToLower(discretization%ctype)
   189)         select case(trim(discretization%ctype))
   190)           case('structured')
   191)             discretization%itype = STRUCTURED_GRID
   192)             call InputReadWord(input,option,structured_grid_ctype,PETSC_TRUE)
   193)             call InputDefaultMsg(input,option,'grid_structured_type') 
   194)             call StringToLower(structured_grid_ctype)
   195)             select case(trim(structured_grid_ctype))
   196)               case('cartesian')
   197)                 structured_grid_itype = CARTESIAN_GRID
   198)               case('cylindrical')
   199)                 structured_grid_itype = CYLINDRICAL_GRID
   200)               case('spherical')
   201)                 structured_grid_itype = SPHERICAL_GRID
   202)               case default
   203)                 structured_grid_itype = CARTESIAN_GRID
   204)                 structured_grid_ctype = 'cartesian'
   205)             end select
   206)           case('unstructured','unstructured_explicit','unstructured_polyhedra')
   207)             discretization%itype = UNSTRUCTURED_GRID
   208)             word = discretization%ctype
   209)             discretization%ctype = 'unstructured'
   210)             select case(word)
   211)               case('unstructured')
   212)                 unstructured_grid_itype = IMPLICIT_UNSTRUCTURED_GRID
   213)                 unstructured_grid_ctype = 'implicit unstructured'
   214)               case('unstructured_explicit')
   215)                 unstructured_grid_itype = EXPLICIT_UNSTRUCTURED_GRID
   216)                 unstructured_grid_ctype = 'explicit unstructured'
   217)               case('unstructured_polyhedra')
   218)                 unstructured_grid_itype = POLYHEDRA_UNSTRUCTURED_GRID
   219)                 unstructured_grid_ctype = 'polyhedra unstructured'
   220)             end select
   221)             call InputReadNChars(input,option,discretization%filename, &
   222)                                  MAXSTRINGLENGTH,PETSC_TRUE)
   223)             call InputErrorMsg(input,option,'unstructured filename','GRID')
   224)           case default
   225)             call InputKeywordUnrecognized(word,'discretization type',option)
   226)         end select    
   227)       case('NXYZ')
   228)         call InputReadInt(input,option,nx)
   229)         call InputErrorMsg(input,option,'nx','GRID')
   230)         call InputReadInt(input,option,ny)
   231)         call InputErrorMsg(input,option,'ny','GRID')
   232)         call InputReadInt(input,option,nz)
   233)         call InputErrorMsg(input,option,'nz','GRID')
   234)         if (structured_grid_itype /= CARTESIAN_GRID) then
   235)           ny = 1 ! cylindrical and spherical have 1 cell in Y
   236)           if (structured_grid_itype /= CYLINDRICAL_GRID) nz = 1 ! spherical has 1 cell in Z
   237)         endif
   238)       case('ORIGIN')
   239)         call InputReadDouble(input,option, &
   240)                              discretization%origin_global(X_DIRECTION))
   241)         call InputErrorMsg(input,option,'X direction','Origin')
   242)         call InputReadDouble(input,option, &
   243)                              discretization%origin_global(Y_DIRECTION))
   244)         call InputErrorMsg(input,option,'Y direction','Origin')
   245)         call InputReadDouble(input,option, &
   246)                              discretization%origin_global(Z_DIRECTION))
   247)         call InputErrorMsg(input,option,'Z direction','Origin')        
   248)       case('FILE','GRAVITY','INVERT_Z','MAX_CELLS_SHARING_A_VERTEX',&
   249)            'STENCIL_WIDTH','STENCIL_TYPE','FLUX_METHOD')
   250)       case('DXYZ','BOUNDS')
   251)         call InputSkipToEND(input,option,word) 
   252)       case default
   253)         call InputKeywordUnrecognized(word,'DISCRETIZATION',option)
   254)     end select 
   255)   enddo  
   256) 
   257)   if (discretization%itype == NULL_GRID) then
   258)     option%io_buffer = 'Discretization type not defined under ' // &
   259)                        'keyword GRID.' 
   260)     call printErrMsg(option)
   261)   endif
   262)   
   263)   grid => GridCreate()
   264)   select case(discretization%itype)
   265)     case(UNSTRUCTURED_GRID)
   266)       un_str_grid => UGridCreate()
   267)       select case(unstructured_grid_itype)
   268)         case(IMPLICIT_UNSTRUCTURED_GRID)
   269)           if (index(discretization%filename,'.h5') > 0) then
   270) #if !defined(PETSC_HAVE_HDF5)
   271)             option%io_buffer = 'PFLOTRAN must be built with HDF5 ' // &
   272)               'support to read unstructured grid .h5 files'
   273)             call printErrMsg(option)
   274) #else
   275) 
   276) #ifdef SCORPIO
   277)             call UGridReadHDF5PIOLib(un_str_grid,discretization%filename,option)
   278) #else
   279)             call UGridReadHDF5(un_str_grid,discretization%filename,option)
   280) #endif
   281) ! #ifdef SCORPIO
   282) 
   283) #endif
   284) !#if !defined(PETSC_HAVE_HDF5)
   285) 
   286)           else
   287)             call UGridRead(un_str_grid,discretization%filename,option)
   288)           endif
   289)           grid%unstructured_grid => un_str_grid
   290)         case(EXPLICIT_UNSTRUCTURED_GRID)
   291)           un_str_grid%explicit_grid => UGridExplicitCreate()
   292)           call UGridExplicitRead(un_str_grid, &
   293)                                  discretization%filename,option)
   294)           grid%unstructured_grid => un_str_grid
   295)         case(POLYHEDRA_UNSTRUCTURED_GRID)
   296)           un_str_grid%polyhedra_grid => UGridPolyhedraCreate()
   297)           if (index(discretization%filename,'.h5') > 0 ) then
   298)             !call UGridPolyhedraReadHDF5(un_str_grid,discretization%filename,option)
   299)             call printErrMsg(option,'Add UGridPolyhedraReadHDF5')
   300)           else
   301)             call UGridPolyhedraRead(un_str_grid,discretization%filename,option)
   302)           endif
   303)           grid%unstructured_grid => un_str_grid
   304)       end select
   305)       grid%itype = unstructured_grid_itype
   306)       grid%ctype = unstructured_grid_ctype
   307)     case(STRUCTURED_GRID)      
   308)       if (nx*ny*nz <= 0) &
   309)         call printErrMsg(option,'NXYZ not set correctly for structured grid.')
   310)       str_grid => StructGridCreate()
   311)       str_grid%nx = nx
   312)       str_grid%ny = ny
   313)       str_grid%nz = nz
   314)       str_grid%nxy = str_grid%nx*str_grid%ny
   315)       str_grid%nmax = str_grid%nxy*str_grid%nz
   316)       grid%structured_grid => str_grid
   317)       grid%nmax = str_grid%nmax
   318)       grid%structured_grid%itype = structured_grid_itype
   319)       grid%structured_grid%ctype = structured_grid_ctype
   320)       grid%itype = discretization%itype
   321)       grid%ctype = discretization%ctype
   322)   end select
   323)   discretization%grid => grid
   324)   nullify(grid)
   325) 
   326) end subroutine DiscretizationReadRequiredCards
   327) 
   328) ! ************************************************************************** !
   329) 
   330) subroutine DiscretizationRead(discretization,input,option)
   331)   ! 
   332)   ! Reads a discretization from the input file
   333)   ! 
   334)   ! Author: Glenn Hammond
   335)   ! Date: 11/01/07
   336)   ! 
   337) 
   338)   use Option_module
   339)   use Input_Aux_module
   340)   use String_module
   341) 
   342)   implicit none
   343) 
   344)   type(option_type), pointer :: option
   345)   type(input_type), pointer :: input
   346)   type(discretization_type),pointer :: discretization
   347)   character(len=MAXWORDLENGTH) :: word
   348)   type(grid_type), pointer :: grid, grid2
   349)   type(grid_structured_type), pointer :: str_grid
   350)   type(grid_unstructured_type), pointer :: un_str_grid
   351)   character(len=MAXWORDLENGTH) :: structured_grid_ctype
   352)   character(len=MAXSTRINGLENGTH) :: filename
   353)   character(len=MAXSTRINGLENGTH) :: string
   354) 
   355)   PetscInt :: structured_grid_itype
   356)   PetscInt :: nx, ny, nz
   357)   PetscInt :: i
   358)   PetscReal :: tempreal
   359) 
   360)   nx = 0
   361)   ny = 0
   362)   nz = 0
   363) 
   364) ! we initialize the word to blanks to avoid error reported by valgrind
   365)   word = ''
   366) 
   367)   do
   368)   
   369)     call InputReadPflotranString(input,option)
   370)     if (input%ierr /= 0) exit
   371) 
   372)     if (InputCheckExit(input,option)) exit
   373) 
   374)     call InputReadWord(input,option,word,PETSC_TRUE)
   375)     call InputErrorMsg(input,option,'keyword','GRID')
   376)     call StringToUpper(word)
   377)       
   378)     select case(trim(word))
   379)       case('TYPE','NXYZ','ORIGIN','FILE')
   380)       case('DXYZ')
   381)         select case(discretization%itype)
   382)           case(STRUCTURED_GRID)
   383)             call StructGridReadDXYZ(discretization%grid%structured_grid,input,option)
   384)           case default
   385)             call printErrMsg(option,&
   386)                            'Keyword "DXYZ" not supported for unstructured grid')
   387)         end select
   388)         call InputReadPflotranString(input,option) ! read END card
   389)         call InputReadStringErrorMsg(input,option,'DISCRETIZATION,DXYZ,END')
   390)         if (.not.(InputCheckExit(input,option))) then
   391)           option%io_buffer = 'Card DXYZ should include either 3 entires ' // &
   392)                    '(one for each grid direction or NX+NY+NZ entries)'
   393)           call printErrMsg(option)
   394)         endif
   395)       case('BOUNDS')
   396)         select case(discretization%itype)
   397)           case(STRUCTURED_GRID)
   398)             grid => discretization%grid
   399) 
   400)             ! read first line and we will split off the legacy approach vs. new
   401)             call InputReadPflotranString(input,option)
   402)             call InputReadStringErrorMsg(input,option, &
   403)                                        'DISCRETIZATION,BOUNDS,Min Coordinates')
   404)             select case(grid%structured_grid%itype)
   405)               case(CARTESIAN_GRID)
   406)                 i = 3
   407)               case(CYLINDRICAL_GRID)
   408)                 i = 2
   409)               case(SPHERICAL_GRID)
   410)                 i = 1
   411)             end select
   412)             call InputReadNDoubles(input,option, &
   413)                                    grid%structured_grid%bounds(:,LOWER), &
   414)                                    i)
   415)             call InputErrorMsg(input,option,'Minimum Coordinate','BOUNDS')
   416)             call InputReadPflotranString(input,option)
   417)             call InputReadStringErrorMsg(input,option, &
   418)                                         'DISCRETIZATION,BOUNDS,Min Coordinates')
   419)             call InputReadNDoubles(input,option, &
   420)                                    grid%structured_grid%bounds(:,UPPER), &
   421)                                    i)
   422)             call InputErrorMsg(input,option,'Maximum Coordinate','BOUNDS')
   423)             if (grid%structured_grid%itype == CYLINDRICAL_GRID) then
   424)               ! 2 values were read in in x and y locations, must move y value
   425)               ! to z as it was really z.
   426)               grid%structured_grid%bounds(Z_DIRECTION,LOWER) = &
   427)                 grid%structured_grid%bounds(Y_DIRECTION,LOWER)
   428)               grid%structured_grid%bounds(Z_DIRECTION,UPPER) = &
   429)                 grid%structured_grid%bounds(Y_DIRECTION,UPPER)
   430)               ! set y bounds to 0 and 1
   431)               grid%structured_grid%bounds(Y_DIRECTION,LOWER) = 0.d0
   432)               grid%structured_grid%bounds(Y_DIRECTION,UPPER) = 1.d0
   433)             endif
   434)             if (grid%structured_grid%itype == SPHERICAL_GRID) then
   435)               grid%structured_grid%bounds(Y_DIRECTION,LOWER) = 0.d0
   436)               grid%structured_grid%bounds(Y_DIRECTION,UPPER) = 1.d0
   437)               grid%structured_grid%bounds(Z_DIRECTION,LOWER) = 0.d0
   438)               grid%structured_grid%bounds(Z_DIRECTION,UPPER) = 1.d0
   439)             endif
   440)             call InputReadPflotranString(input,option)
   441)             call InputReadStringErrorMsg(input,option, &
   442)                                          'DISCRETIZATION,BOUNDS,END')
   443)             if (.not.(InputCheckExit(input,option))) then
   444)               if (OptionPrintToScreen(option)) then
   445)                 if (grid%structured_grid%itype == CARTESIAN_GRID) then
   446)                   print *, 'BOUNDS card for a cartesian structured grid ' // &
   447)                     'must include 4 lines.  I.e.'
   448)                   print *, 'BOUNDS'
   449)                   print *, '  x_min  y_min  z_min'
   450)                   print *, '  x_max  y_max  z_max'
   451)                   print *, 'END'
   452)                 else if (grid%structured_grid%itype == CYLINDRICAL_GRID) then
   453)                   print *, 'BOUNDS card for a cylindrical structured grid ' // &
   454)                     'must include 4 lines.  I.e.'
   455)                   print *, 'BOUNDS'
   456)                   print *, '  r_min  z_min'
   457)                   print *, '  r_max  z_max'
   458)                   print *, 'END'
   459)                 else if (grid%structured_grid%itype == SPHERICAL_GRID) then
   460)                   print *, 'BOUNDS card for a spherical structured grid ' // &
   461)                     'must include 4 lines.  I.e.'
   462)                   print *, 'BOUNDS'
   463)                   print *, '  r_min'
   464)                   print *, '  r_max'
   465)                   print *, 'END'
   466)                 endif
   467)               endif
   468)               stop
   469)             endif            
   470)             discretization%origin_global(X_DIRECTION) = &
   471)               grid%structured_grid%bounds(X_DIRECTION,LOWER)
   472)             discretization%origin_global(Y_DIRECTION) = &
   473)               grid%structured_grid%bounds(Y_DIRECTION,LOWER)
   474)             discretization%origin_global(Z_DIRECTION) = &
   475)               grid%structured_grid%bounds(Z_DIRECTION,LOWER)
   476)         end select
   477)       case ('GRAVITY')
   478)         call InputReadDouble(input,option,option%gravity(X_DIRECTION))
   479)         call InputErrorMsg(input,option,'x-direction','GRAVITY')
   480)         call InputReadDouble(input,option,option%gravity(Y_DIRECTION))
   481)         call InputErrorMsg(input,option,'y-direction','GRAVITY')
   482)         call InputReadDouble(input,option,option%gravity(Z_DIRECTION))
   483)         call InputErrorMsg(input,option,'z-direction','GRAVITY')
   484)         if (option%myrank == option%io_rank .and. &
   485)             option%print_to_screen) &
   486)           write(option%fid_out,'(/," *GRAV",/, &
   487)             & "  gravity    = "," [m/s^2]",3x,1p3e12.4 &
   488)             & )') option%gravity(1:3)
   489)       case ('MAX_CELLS_SHARING_A_VERTEX')
   490)         if (associated(discretization%grid%unstructured_grid)) then
   491)           call InputReadInt(input,option,discretization%grid% &
   492)                             unstructured_grid%max_cells_sharing_a_vertex)
   493)           call InputErrorMsg(input,option,'max_cells_sharing_a_vertex', &
   494)                              'GRID')
   495)         endif          
   496)       case ('INVERT_Z')
   497)       case ('STENCIL_WIDTH')
   498)         call InputReadInt(input,option,discretization%stencil_width)
   499)         call InputErrorMsg(input,option,'stencil_width', &
   500)                            'GRID')
   501)       case ('STENCIL_TYPE')
   502)         call InputReadWord(input,option,word,PETSC_TRUE)
   503)         call InputErrorMsg(input,option,'keyword','GRID')
   504)         call StringToUpper(word)
   505)         select case(trim(word))
   506)           case ('BOX')
   507)             discretization%stencil_type = DMDA_STENCIL_BOX
   508)           case ('STAR')
   509)             discretization%stencil_type = DMDA_STENCIL_STAR
   510)           case default
   511)             call InputKeywordUnrecognized(word, &
   512)                    'DISCRETIZATION,stencil type',option)
   513)         end select
   514)       case default
   515)         call InputKeywordUnrecognized(word,'GRID',option)
   516)     end select 
   517)   enddo  
   518) 
   519)   select case(discretization%itype)
   520)     case(STRUCTURED_GRID)
   521)       if (discretization%grid%structured_grid%invert_z_axis) then
   522)         option%gravity(Z_DIRECTION) = -1.d0*option%gravity(Z_DIRECTION)
   523)       endif
   524)   end select
   525)   
   526) end subroutine DiscretizationRead
   527) 
   528) ! ************************************************************************** !
   529) 
   530) subroutine DiscretizationCreateDMs(discretization, o_nflowdof, o_ntrandof, &
   531)                                     o_nphase, o_ngeomechdof, o_n_stress_strain_dof, option)
   532) 
   533)   ! 
   534)   ! creates distributed, parallel meshes/grids
   535)   ! If there are multiple degrees of freedom per grid cell, this will call
   536)   ! DiscretizationCreateDM() multiple times to create the DMs corresponding
   537)   ! to one degree of freedom grid cell and those corresponding to multiple
   538)   ! degrees of freedom per cell.
   539)   ! 
   540)   ! Author: Glenn Hammond
   541)   ! Date: 02/08/08
   542)   ! 
   543)       
   544)   use Option_module    
   545)       
   546)   implicit none
   547)   
   548)   type(discretization_type) :: discretization
   549)   PetscInt, intent(in) :: o_nflowdof
   550)   PetscInt, intent(in) :: o_ntrandof
   551)   PetscInt, intent(in) :: o_nphase
   552)   PetscInt, intent(in) :: o_ngeomechdof
   553)   PetscInt, intent(in) :: o_n_stress_strain_dof
   554)   type(option_type) :: option
   555)       
   556)   PetscInt :: ndof
   557)   !PetscInt, parameter :: stencil_width = 1
   558)   PetscErrorCode :: ierr
   559)   PetscInt :: i
   560)   type(grid_unstructured_type), pointer :: ugrid
   561) 
   562)   select case(discretization%itype)
   563)     case(STRUCTURED_GRID)
   564)       discretization%dm_index_to_ndof(ONEDOF) = 1
   565)       discretization%dm_index_to_ndof(NPHASEDOF) = o_nphase
   566)       discretization%dm_index_to_ndof(NFLOWDOF) = o_nflowdof
   567)       discretization%dm_index_to_ndof(NTRANDOF) = o_ntrandof
   568)     case(UNSTRUCTURED_GRID)
   569) 
   570)     
   571)       select case(discretization%grid%itype)
   572)         case(IMPLICIT_UNSTRUCTURED_GRID)
   573)           ! petsc will call parmetis to calculate the graph/dual
   574) #if !defined(PETSC_HAVE_PARMETIS)
   575)             option%io_buffer = &
   576)              'Must compile with Parmetis in order to use implicit unstructured grids.'
   577)             call printErrMsg(option)
   578) #endif
   579)           call UGridDecompose(discretization%grid%unstructured_grid, &
   580)                               option)
   581)         case(EXPLICIT_UNSTRUCTURED_GRID)
   582) #if !defined(PETSC_HAVE_PARMETIS) && !defined(PETSC_HAVE_PTSCOTCH)
   583)             option%io_buffer = &
   584)              'Must compile with either Parmetis or PTSCOTCH in order to use explicit unstructured grids.'
   585)             call printErrMsg(option)
   586) #endif
   587)           ugrid => discretization%grid%unstructured_grid
   588)           call UGridExplicitDecompose(ugrid,option)
   589)         case(POLYHEDRA_UNSTRUCTURED_GRID)
   590) #if !defined(PETSC_HAVE_PARMETIS)
   591)             option%io_buffer = &
   592)              'Must compile with Parmetis in order to use implicit unstructured grids.'
   593)             call printErrMsg(option)
   594) #endif
   595)           ugrid => discretization%grid%unstructured_grid
   596)           call UGridPolyhedraDecompose(ugrid,option)
   597)       end select
   598)   end select
   599) 
   600) 
   601)   !-----------------------------------------------------------------------
   602)   ! Generate the DM objects that will manage communication.
   603)   !-----------------------------------------------------------------------
   604)   ndof = 1
   605)   call DiscretizationCreateDM(discretization,discretization%dm_1dof, &
   606)                               ndof,discretization%stencil_width, &
   607)                               discretization%stencil_type,option)
   608)   
   609)   if (o_nflowdof > 0) then
   610)     ndof = o_nflowdof
   611)     call DiscretizationCreateDM(discretization,discretization%dm_nflowdof, &
   612)                                 ndof,discretization%stencil_width, &
   613)                                 discretization%stencil_type,option)
   614)   endif
   615)   
   616)   if (o_ntrandof > 0) then
   617)     ndof = o_ntrandof
   618)     call DiscretizationCreateDM(discretization,discretization%dm_ntrandof, &
   619)                                 ndof,discretization%stencil_width, &
   620)                                 discretization%stencil_type,option)
   621)   endif
   622) 
   623)   if (o_ngeomechdof > 0) then
   624)     ndof = o_n_stress_strain_dof
   625)     call DiscretizationCreateDM(discretization,discretization%dm_n_stress_strain_dof, &
   626)                                 ndof,discretization%stencil_width, &
   627)                                 discretization%stencil_type,option)
   628)   endif
   629) 
   630)   select case(discretization%itype)
   631)     case(STRUCTURED_GRID)
   632)       ! this function must be called to set up str_grid%lxs, etc.
   633)       call StructGridComputeLocalBounds(discretization%grid%structured_grid, &
   634)                                         discretization%dm_1dof%dm,option)    
   635)       discretization%grid%nlmax = discretization%grid%structured_grid%nlmax
   636)       discretization%grid%ngmax = discretization%grid%structured_grid%ngmax
   637)       discretization%grid%global_offset = &
   638)         discretization%grid%structured_grid%global_offset
   639)     case(UNSTRUCTURED_GRID)
   640)       discretization%grid%nmax = discretization%grid%unstructured_grid%nmax
   641)       discretization%grid%nlmax = discretization%grid%unstructured_grid%nlmax
   642)       discretization%grid%ngmax = discretization%grid%unstructured_grid%ngmax
   643)       discretization%grid%global_offset = &
   644)         discretization%grid%unstructured_grid%global_offset
   645)   end select
   646) 
   647) end subroutine DiscretizationCreateDMs
   648) 
   649) ! ************************************************************************** !
   650) 
   651) subroutine DiscretizationCreateDM(discretization,dm_ptr,ndof,stencil_width, &
   652)                                   stencil_type,option)
   653)   ! 
   654)   ! creates a distributed, parallel mesh/grid
   655)   ! 
   656)   ! Author: Glenn Hammond
   657)   ! Date: 02/08/08
   658)   ! 
   659) 
   660)   use Option_module
   661)   
   662)   implicit none
   663)   
   664)   type(discretization_type) :: discretization
   665)   type(dm_ptr_type), pointer :: dm_ptr
   666)   PetscInt :: ndof
   667)   PetscInt :: stencil_width
   668)   PetscEnum :: stencil_type
   669)   type(option_type) :: option
   670)   PetscErrorCode :: ierr
   671) 
   672)   select case(discretization%itype)
   673)     case(STRUCTURED_GRID)
   674)       call StructGridCreateDM(discretization%grid%structured_grid, &
   675)                               dm_ptr%dm,ndof,stencil_width,stencil_type,option)
   676)     case(UNSTRUCTURED_GRID)
   677)       call UGridCreateUGDMShell(discretization%grid%unstructured_grid, &
   678)                            dm_ptr%dm,dm_ptr%ugdm,ndof,option)
   679)   end select
   680) 
   681) end subroutine DiscretizationCreateDM
   682) 
   683) ! ************************************************************************** !
   684) 
   685) subroutine DiscretizationCreateVector(discretization,dm_index,vector, &
   686)                                       vector_type,option)
   687)   ! 
   688)   ! Creates a global PETSc vector
   689)   ! 
   690)   ! Author: Glenn Hammond
   691)   ! Date: 10/24/07
   692)   ! 
   693)   use Option_module                                      
   694) 
   695)   implicit none
   696)   
   697)   type(discretization_type) :: discretization
   698)   PetscInt :: dm_index
   699)   Vec :: vector
   700)   PetscInt :: vector_type
   701)   type(option_type) :: option
   702)   PetscInt :: ndof
   703)   PetscErrorCode :: ierr
   704)   
   705)   type(dm_ptr_type), pointer :: dm_ptr
   706)   
   707)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
   708) 
   709)   select case (vector_type)
   710)     case(GLOBAL)
   711)       call DMCreateGlobalVector(dm_ptr%dm,vector,ierr);CHKERRQ(ierr)
   712)     case(LOCAL)
   713)       call DMCreateLocalVector(dm_ptr%dm,vector,ierr);CHKERRQ(ierr)
   714)     case(NATURAL)
   715)       select case(discretization%itype)
   716)         case(STRUCTURED_GRID)
   717)           call DMDACreateNaturalVector(dm_ptr%dm,vector,ierr);CHKERRQ(ierr)
   718)         case(UNSTRUCTURED_GRID)
   719)           call UGridDMCreateVector(discretization%grid%unstructured_grid, &
   720)                                    dm_ptr%ugdm,vector, &
   721)                                    vector_type,option)
   722)         end select
   723)   end select
   724) 
   725)   call VecSet(vector,0.d0,ierr);CHKERRQ(ierr)
   726)   
   727) end subroutine DiscretizationCreateVector
   728) 
   729) ! ************************************************************************** !
   730) 
   731) subroutine DiscretizationDuplicateVector(discretization,vector1,vector2)
   732)   ! 
   733)   ! Creates a global PETSc vector
   734)   ! 
   735)   ! Author: Glenn Hammond
   736)   ! Date: 10/24/07
   737)   ! 
   738) 
   739)   implicit none
   740)   
   741)   type(discretization_type) :: discretization
   742)   Vec :: vector1
   743)   Vec :: vector2
   744)   
   745)   PetscErrorCode :: ierr
   746)   call VecDuplicate(vector1,vector2,ierr);CHKERRQ(ierr)
   747)   call VecCopy(vector1,vector2,ierr);CHKERRQ(ierr)
   748)   
   749) end subroutine DiscretizationDuplicateVector
   750) 
   751) ! ************************************************************************** !
   752) 
   753) function DiscretizationGetDMPtrFromIndex(discretization,dm_index)
   754)   ! 
   755)   ! Returns the integer pointer for the DM referenced
   756)   ! 
   757)   ! Author: Glenn Hammond
   758)   ! Date: 02/08/08
   759)   ! 
   760) 
   761)   implicit none
   762)   
   763)   type(discretization_type) :: discretization
   764)   PetscInt :: dm_index
   765)   
   766)   type(dm_ptr_type), pointer :: DiscretizationGetDMPtrFromIndex
   767)   
   768)   select case (dm_index)
   769)     case(ONEDOF)
   770)       DiscretizationGetDMPtrFromIndex => discretization%dm_1dof
   771)     case(NFLOWDOF)
   772)       DiscretizationGetDMPtrFromIndex => discretization%dm_nflowdof
   773)     case(NTRANDOF)
   774)       DiscretizationGetDMPtrFromIndex => discretization%dm_ntrandof
   775)     case(NGEODOF)
   776)       DiscretizationGetDMPtrFromIndex => discretization%dm_n_stress_strain_dof
   777)   end select  
   778)   
   779) end function DiscretizationGetDMPtrFromIndex
   780) 
   781) ! ************************************************************************** !
   782) 
   783) function DiscretizationGetDMCPtrFromIndex(discretization,dm_index)
   784) 
   785)   implicit none
   786)   
   787)   type(discretization_type) :: discretization
   788)   PetscInt :: dm_index
   789)   
   790)   type(dm_ptr_type), pointer :: DiscretizationGetDMCPtrFromIndex(:)
   791)   
   792)   select case (dm_index)
   793)     case(NFLOWDOF)
   794)       DiscretizationGetDMCPtrFromIndex => discretization%dmc_nflowdof
   795)     case(NTRANDOF)
   796)       DiscretizationGetDMCPtrFromIndex => discretization%dmc_ntrandof
   797)   end select  
   798)   
   799) end function DiscretizationGetDMCPtrFromIndex
   800) 
   801) ! ************************************************************************** !
   802) 
   803) subroutine DiscretizationCreateJacobian(discretization,dm_index,mat_type,Jacobian,option)
   804)   ! 
   805)   ! Creates Jacobian matrix associated with discretization
   806)   ! 
   807)   ! Author: Glenn Hammond
   808)   ! Date: 10/24/07
   809)   ! 
   810) 
   811)   use Option_module
   812)   
   813)   implicit none
   814)   
   815) #include "petsc/finclude/petscis.h"
   816) #include "petsc/finclude/petscis.h90"
   817) 
   818)   type(discretization_type) :: discretization
   819)   PetscInt :: dm_index
   820)   PetscErrorCode :: ierr
   821)   MatType :: mat_type
   822)   Mat :: Jacobian
   823)   type(option_type) :: option
   824)   PetscInt :: ndof, stencilsize
   825)   PetscInt, pointer :: indices(:)
   826)   PetscInt :: ngmax
   827)   PetscInt :: imax, nlevels, ln, npatches, pn, i
   828)   type(dm_ptr_type), pointer :: dm_ptr
   829)   ISLocalToGlobalMapping :: ptmap
   830)   PetscInt :: islocal
   831) 
   832)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
   833) 
   834) 
   835)   select case(discretization%itype)
   836)     case(STRUCTURED_GRID)
   837)       call DMSetMatType(dm_ptr%dm,mat_type,ierr);CHKERRQ(ierr)
   838)       call DMCreateMatrix(dm_ptr%dm,Jacobian,ierr);CHKERRQ(ierr)
   839)     case(UNSTRUCTURED_GRID)
   840)       call UGridDMCreateJacobian(discretization%grid%unstructured_grid, &
   841)                                  dm_ptr%ugdm,mat_type,Jacobian,option)
   842)   end select
   843)   call MatSetOption(Jacobian,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE, &
   844)                     ierr);CHKERRQ(ierr)
   845)   call MatSetOption(Jacobian,MAT_ROW_ORIENTED,PETSC_FALSE,ierr);CHKERRQ(ierr)
   846)   call MatSetOption(Jacobian,MAT_NO_OFF_PROC_ZERO_ROWS,PETSC_TRUE, &
   847)                     ierr);CHKERRQ(ierr)
   848) 
   849) end subroutine DiscretizationCreateJacobian
   850) 
   851) ! ************************************************************************** !
   852) 
   853) subroutine DiscretizationCreateInterpolation(discretization,dm_index, &
   854)                                              interpolation,mg_levels_x, &
   855)                                              mg_levels_y, mg_levels_z, &
   856)                                              option)
   857)   ! 
   858)   ! Creates interpolation matrix associated
   859)   ! with discretization for geometric multigrid.
   860)   ! 
   861)   ! Author: Richard Mills
   862)   ! Date: 4/25/08.
   863)   ! 
   864) 
   865)   use Option_module
   866)   
   867)   implicit none
   868)   
   869)   type(discretization_type) :: discretization
   870)   PetscInt :: dm_index
   871)   PetscErrorCode :: ierr
   872)   Mat, pointer :: interpolation(:)
   873)   PetscInt :: mg_levels_x, mg_levels_y, mg_levels_z
   874)   type(option_type) :: option
   875) 
   876)   PetscInt :: mg_levels
   877)   PetscInt :: refine_x, refine_y, refine_z
   878)   type(dm_ptr_type), pointer :: dm_ptr
   879)   type(dm_ptr_type), pointer :: dmc_ptr(:)
   880)   PetscInt :: i
   881)   type(dm_ptr_type), pointer :: dm_fine_ptr
   882)     ! Used to point to finer-grid DM in the loop that constructst the 
   883)     ! interpolation hierarchy.
   884)   
   885)   mg_levels = max(mg_levels_x, mg_levels_y, mg_levels_z)
   886)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
   887) !  dmc_ptr = DiscretizationGetDMCPtrFromIndex(discretization,dm_index)
   888)   select case (dm_index)
   889)     case(NFLOWDOF)
   890)       allocate(discretization%dmc_nflowdof(mg_levels))
   891)       do i=1, mg_levels
   892)         discretization%dmc_nflowdof(i)%dm = 0
   893)         nullify(discretization%dmc_nflowdof(i)%ugdm)
   894)       enddo
   895)       dmc_ptr => discretization%dmc_nflowdof
   896)     case(NTRANDOF)
   897)       allocate(discretization%dmc_ntrandof(mg_levels))
   898)       do i=1, mg_levels
   899)         discretization%dmc_ntrandof(i)%dm = 0
   900)         nullify(discretization%dmc_ntrandof(i)%ugdm)
   901)       enddo
   902)       dmc_ptr => discretization%dmc_ntrandof
   903)   end select  
   904)    
   905)   allocate(interpolation(mg_levels))
   906) 
   907)   select case(discretization%itype)
   908)     case(STRUCTURED_GRID)
   909)       dm_fine_ptr => dm_ptr
   910)       refine_x = 2; refine_y = 2; refine_z = 2
   911)       do i=mg_levels-1,1,-1
   912)         ! If number of coarsenings performed so far exceeds mg_levels_x-1, 
   913)         ! set refine_x = 1; likewise for y and z.
   914)         if (i <= mg_levels - mg_levels_x ) refine_x = 1
   915)         if (i <= mg_levels - mg_levels_y ) refine_y = 1
   916)         if (i <= mg_levels - mg_levels_z ) refine_z = 1
   917)         call DMDASetRefinementFactor(dm_fine_ptr%dm, refine_x, refine_y, refine_z, &
   918)                                    ierr);CHKERRQ(ierr)
   919)         call DMDASetInterpolationType(dm_fine_ptr%dm, DMDA_Q0,  &
   920)                                       ierr);CHKERRQ(ierr)
   921)         call DMCoarsen(dm_fine_ptr%dm, option%mycomm, dmc_ptr(i)%dm,  &
   922)                        ierr);CHKERRQ(ierr)
   923) #ifndef DMGET
   924)         call DMCreateInterpolation(dmc_ptr(i)%dm, dm_fine_ptr%dm, &
   925)                                    interpolation(i), PETSC_NULL_OBJECT,  &
   926)                                    ierr);CHKERRQ(ierr)
   927) #else
   928)         call DMGetInterpolation(dmc_ptr(i)%dm, dm_fine_ptr%dm, &
   929)                                 interpolation(i), PETSC_NULL_OBJECT,  &
   930)                                 ierr);CHKERRQ(ierr)
   931) #endif
   932)         dm_fine_ptr => dmc_ptr(i)
   933)       enddo
   934)     case(UNSTRUCTURED_GRID)
   935)   end select
   936) 
   937) end subroutine DiscretizationCreateInterpolation
   938) 
   939) ! ************************************************************************** !
   940) 
   941) subroutine DiscretizationCreateColoring(discretization,dm_index,option,coloring)
   942)   ! 
   943)   ! Creates ISColoring for discretization
   944)   ! 
   945)   ! Author: Glenn Hammond
   946)   ! Date: 10/24/07
   947)   ! 
   948) 
   949)   use Option_module
   950)   
   951)   implicit none
   952) 
   953) #include "petsc/finclude/petscis.h"
   954) #include "petsc/finclude/petscis.h90"
   955)   
   956)   type(discretization_type) :: discretization
   957)   PetscInt :: dm_index
   958)   PetscErrorCode :: ierr
   959)   type(option_type) :: option
   960)   ISColoring :: coloring
   961) 
   962)   type(dm_ptr_type), pointer :: dm_ptr
   963)   
   964)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
   965)     
   966)   select case(discretization%itype)
   967)     case(STRUCTURED_GRID)
   968) #ifndef DMGET
   969)       call DMCreateColoring(dm_ptr%dm,IS_COLORING_GLOBAL,MATBAIJ,coloring,&
   970)                             ierr);CHKERRQ(ierr)
   971) #else
   972)       call DMGetColoring(dm_ptr%dm,IS_COLORING_GLOBAL,MATBAIJ,coloring, &
   973)                          ierr);CHKERRQ(ierr)
   974) #endif
   975)       ! I have set the above to use matrix type MATBAIJ, as that is what we 
   976)       ! usually want (note: for DAs with 1 degree of freedom per grid cell, 
   977)       ! the MATAIJ and MATBAIJ colorings should be equivalent).  What we should 
   978)       ! eventually do here is query the type of the Jacobian matrix, but I'm 
   979)       ! not sure of the best way to do this, as this is currently stashed in 
   980)       ! the 'solver' object. --RTM
   981)     case(UNSTRUCTURED_GRID)
   982)   end select
   983)   
   984) end subroutine DiscretizationCreateColoring
   985) 
   986) ! ************************************************************************** !
   987) 
   988) subroutine DiscretizationGlobalToLocal(discretization,global_vec,local_vec,dm_index)
   989)   ! 
   990)   ! Performs global to local communication with DM
   991)   ! Note that 'dm_index' should correspond to one of the macros defined
   992)   ! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
   993)   ! 
   994)   ! Author: Glenn Hammond
   995)   ! Date: 10/24/07
   996)   ! 
   997) 
   998)   implicit none
   999) 
  1000)   type(discretization_type) :: discretization
  1001)   Vec :: global_vec
  1002)   Vec :: local_vec
  1003)   PetscInt :: dm_index
  1004)   PetscErrorCode :: ierr
  1005)   type(dm_ptr_type), pointer :: dm_ptr
  1006)   
  1007)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1008)     
  1009)   call DMGlobalToLocalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec, &
  1010)                             ierr);CHKERRQ(ierr)
  1011)   call DMGlobalToLocalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec, &
  1012)                           ierr);CHKERRQ(ierr)
  1013)   
  1014) end subroutine DiscretizationGlobalToLocal
  1015) 
  1016) ! ************************************************************************** !
  1017) 
  1018) subroutine DiscretizationLocalToGlobal(discretization,local_vec,global_vec,dm_index)
  1019)   ! 
  1020)   ! Performs local to global communication with DM
  1021)   ! Note that 'dm_index' should correspond to one of the macros defined
  1022)   ! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
  1023)   ! 
  1024)   ! Author: Glenn Hammond
  1025)   ! Date: 1/02/08
  1026)   ! 
  1027) 
  1028)   implicit none
  1029)   
  1030)   type(discretization_type) :: discretization
  1031)   Vec :: local_vec
  1032)   Vec :: global_vec
  1033)   PetscInt :: dm_index
  1034)   PetscErrorCode :: ierr
  1035)   type(dm_ptr_type), pointer :: dm_ptr
  1036)   
  1037)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1038)   
  1039)   call DMLocalToGlobalBegin(dm_ptr%dm,local_vec,INSERT_VALUES,global_vec, &
  1040)                             ierr);CHKERRQ(ierr)
  1041)   call DMLocalToGlobalEnd(dm_ptr%dm,local_vec,INSERT_VALUES,global_vec, &
  1042)                           ierr);CHKERRQ(ierr)
  1043)  
  1044) end subroutine DiscretizationLocalToGlobal
  1045) 
  1046) ! ************************************************************************** !
  1047) 
  1048) subroutine DiscretizationLocalToLocal(discretization,local_vec1,local_vec2,dm_index)
  1049)   ! 
  1050)   ! Performs local to local communication with DM
  1051)   ! Some clarification:
  1052)   ! A "local to local" operation, in PETSc parlance, refers to communicating
  1053)   ! values from a local ghosted vector (in which the ghost points are
  1054)   ! irrelevant) and putting those values directly into another ghosted local
  1055)   ! vector (in which those ghost points are set correctly).
  1056)   ! This uses the same communication pattern as a "global to local" operation,
  1057)   ! but a in a "global to local", the originating vector is a PETSc global
  1058)   ! vector, not a ghosted local vector.
  1059)   ! 
  1060)   ! Author: Glenn Hammond
  1061)   ! Date: 11/14/07
  1062)   ! 
  1063) 
  1064)   implicit none
  1065)   
  1066)   type(discretization_type) :: discretization
  1067)   Vec :: local_vec1
  1068)   Vec :: local_vec2
  1069)   PetscInt :: dm_index
  1070)   PetscErrorCode :: ierr
  1071)   type(dm_ptr_type), pointer :: dm_ptr
  1072)   
  1073)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1074)   
  1075)   call DMLocalToLocalBegin(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2, &
  1076)                            ierr);CHKERRQ(ierr)
  1077)   call DMLocalToLocalEnd(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2, &
  1078)                          ierr);CHKERRQ(ierr)
  1079)   
  1080) end subroutine DiscretizationLocalToLocal
  1081) 
  1082) ! ************************************************************************** !
  1083) 
  1084) subroutine DiscretizationGlobalToNatural(discretization,global_vec,natural_vec,dm_index)
  1085)   ! 
  1086)   ! Performs global to natural communication with DM
  1087)   ! 
  1088)   ! Author: Glenn Hammond
  1089)   ! Date: 10/24/07
  1090)   ! 
  1091) 
  1092)   implicit none
  1093)   
  1094)   type(discretization_type) :: discretization
  1095)   Vec :: global_vec
  1096)   Vec :: natural_vec
  1097)   PetscInt :: dm_index
  1098)   PetscErrorCode :: ierr
  1099)   type(dm_ptr_type), pointer :: dm_ptr
  1100)   
  1101)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1102)   
  1103)   select case(discretization%itype)
  1104)     case(STRUCTURED_GRID)
  1105)       call DMDAGlobalToNaturalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec, &
  1106)                                     ierr);CHKERRQ(ierr)
  1107)       call DMDAGlobalToNaturalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec, &
  1108)                                   ierr);CHKERRQ(ierr)
  1109)     case(UNSTRUCTURED_GRID)
  1110)       call VecScatterBegin(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
  1111)                            INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
  1112)       call VecScatterEnd(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
  1113)                          INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
  1114)   end select
  1115)   
  1116) end subroutine DiscretizationGlobalToNatural
  1117) 
  1118) ! ************************************************************************** !
  1119) 
  1120) subroutine DiscretizationNaturalToGlobal(discretization,natural_vec,global_vec,dm_index)
  1121)   ! 
  1122)   ! Performs natural to global communication with DM
  1123)   ! 
  1124)   ! Author: Glenn Hammond
  1125)   ! Date: 01/12/08
  1126)   ! 
  1127) 
  1128)   implicit none
  1129)   
  1130)   type(discretization_type) :: discretization
  1131)   Vec :: natural_vec
  1132)   Vec :: global_vec
  1133)   PetscInt :: dm_index
  1134)   PetscErrorCode :: ierr
  1135)   type(dm_ptr_type), pointer :: dm_ptr
  1136)   
  1137)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1138)   
  1139)   select case(discretization%itype)
  1140)     case(STRUCTURED_GRID)
  1141)       call DMDANaturalToGlobalBegin(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec, &
  1142)                                     ierr);CHKERRQ(ierr)
  1143)       call DMDANaturalToGlobalEnd(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec, &
  1144)                                   ierr);CHKERRQ(ierr)
  1145)     case(UNSTRUCTURED_GRID)
  1146)       call VecScatterBegin(dm_ptr%ugdm%scatter_gton,natural_vec,global_vec, &
  1147)                            INSERT_VALUES,SCATTER_REVERSE,ierr);CHKERRQ(ierr)
  1148)       call VecScatterEnd(dm_ptr%ugdm%scatter_gton,natural_vec,global_vec, &
  1149)                          INSERT_VALUES,SCATTER_REVERSE,ierr);CHKERRQ(ierr)
  1150)   end select
  1151)   
  1152) end subroutine DiscretizationNaturalToGlobal
  1153) 
  1154) ! ************************************************************************** !
  1155) 
  1156) subroutine DiscretizationGlobalToLocalBegin(discretization,global_vec,local_vec,dm_index)
  1157)   ! 
  1158)   ! Begins global to local communication with DM
  1159)   ! Note that 'dm_index' should correspond to one of the macros defined
  1160)   ! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
  1161)   ! 
  1162)   ! Author: Glenn Hammond
  1163)   ! Date: 10/24/07
  1164)   ! 
  1165) 
  1166) 
  1167) 
  1168)   implicit none
  1169)   
  1170)   type(discretization_type) :: discretization
  1171)   Vec :: global_vec
  1172)   Vec :: local_vec
  1173)   PetscInt :: dm_index
  1174)   PetscErrorCode :: ierr
  1175)   type(dm_ptr_type), pointer :: dm_ptr
  1176)   
  1177)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1178)   
  1179)   call DMGlobalToLocalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec, &
  1180)                             ierr);CHKERRQ(ierr)
  1181)   
  1182) end subroutine DiscretizationGlobalToLocalBegin
  1183) 
  1184) ! ************************************************************************** !
  1185) 
  1186) subroutine DiscretizationGlobalToLocalEnd(discretization,global_vec,local_vec,dm_index)
  1187)   ! 
  1188)   ! Ends global to local communication with DM
  1189)   ! Note that 'dm_index' should correspond to one of the macros defined
  1190)   ! in 'definitions.h' such as ONEDOF, NPHASEDOF, etc.  --RTM
  1191)   ! 
  1192)   ! Author: Glenn Hammond
  1193)   ! Date: 10/24/07
  1194)   ! 
  1195) 
  1196)  
  1197) 
  1198)  implicit none
  1199)   
  1200)   type(discretization_type) :: discretization
  1201)   Vec :: global_vec
  1202)   Vec :: local_vec
  1203)   PetscInt :: dm_index
  1204)   PetscErrorCode :: ierr
  1205)   type(dm_ptr_type), pointer :: dm_ptr
  1206)   
  1207)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1208)   
  1209)   call DMGlobalToLocalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,local_vec, &
  1210)                           ierr);CHKERRQ(ierr)
  1211)  
  1212) end subroutine DiscretizationGlobalToLocalEnd
  1213) 
  1214) ! ************************************************************************** !
  1215) 
  1216) subroutine DiscretizationLocalToLocalBegin(discretization,local_vec1,local_vec2,dm_index)
  1217)   ! 
  1218)   ! Begins local to local communication with DM
  1219)   ! 
  1220)   ! Author: Glenn Hammond
  1221)   ! Date: 11/14/07
  1222)   ! 
  1223) 
  1224)   implicit none
  1225)   
  1226)   type(discretization_type) :: discretization
  1227)   Vec :: local_vec1
  1228)   Vec :: local_vec2
  1229)   PetscInt :: dm_index
  1230)   PetscErrorCode :: ierr
  1231)   type(dm_ptr_type), pointer :: dm_ptr
  1232)   
  1233)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1234)   
  1235)   call DMLocalToLocalBegin(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2, &
  1236)                            ierr);CHKERRQ(ierr)
  1237) 
  1238) end subroutine DiscretizationLocalToLocalBegin
  1239) 
  1240) ! ************************************************************************** !
  1241) 
  1242) subroutine DiscretizationLocalToLocalEnd(discretization,local_vec1,local_vec2,dm_index)
  1243)   ! 
  1244)   ! Ends local to local communication with DM
  1245)   ! 
  1246)   ! Author: Glenn Hammond
  1247)   ! Date: 11/14/07
  1248)   ! 
  1249) 
  1250)   implicit none
  1251)   
  1252)   type(discretization_type) :: discretization
  1253)   Vec :: local_vec1
  1254)   Vec :: local_vec2
  1255)   PetscInt :: dm_index
  1256)   PetscErrorCode :: ierr
  1257)   type(dm_ptr_type), pointer :: dm_ptr
  1258)   
  1259)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1260)   
  1261)   call DMLocalToLocalEnd(dm_ptr%dm,local_vec1,INSERT_VALUES,local_vec2, &
  1262)                          ierr);CHKERRQ(ierr)
  1263) 
  1264) end subroutine DiscretizationLocalToLocalEnd
  1265) 
  1266) ! ************************************************************************** !
  1267) 
  1268) subroutine DiscretizGlobalToNaturalBegin(discretization,global_vec,natural_vec,dm_index)
  1269)   ! 
  1270)   ! Begins global to natural communication with DM
  1271)   ! 
  1272)   ! Author: Glenn Hammond
  1273)   ! Date: 10/24/07
  1274)   ! 
  1275) 
  1276)   implicit none
  1277)   
  1278)   type(discretization_type) :: discretization
  1279)   Vec :: global_vec
  1280)   Vec :: natural_vec
  1281)   PetscInt :: dm_index
  1282)   PetscErrorCode :: ierr
  1283)   type(dm_ptr_type), pointer :: dm_ptr
  1284)   
  1285)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1286)   
  1287)   select case(discretization%itype)
  1288)     case(STRUCTURED_GRID)
  1289)       call DMDAGlobalToNaturalBegin(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec, &
  1290)                                     ierr);CHKERRQ(ierr)
  1291)     case(UNSTRUCTURED_GRID)
  1292)       call VecScatterBegin(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
  1293)                            INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
  1294)   end select
  1295)   
  1296) end subroutine DiscretizGlobalToNaturalBegin
  1297) 
  1298) ! ************************************************************************** !
  1299) 
  1300) subroutine DiscretizGlobalToNaturalEnd(discretization,global_vec,natural_vec,dm_index)
  1301)   ! 
  1302)   ! Ends global to natural communication with DM
  1303)   ! 
  1304)   ! Author: Glenn Hammond
  1305)   ! Date: 10/24/07
  1306)   ! 
  1307) 
  1308)   implicit none
  1309)   
  1310)   type(discretization_type) :: discretization
  1311)   Vec :: global_vec
  1312)   Vec :: natural_vec
  1313)   PetscInt :: dm_index
  1314)   PetscErrorCode :: ierr
  1315)   type(dm_ptr_type), pointer :: dm_ptr
  1316)   
  1317)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1318)   
  1319)   select case(discretization%itype)
  1320)     case(STRUCTURED_GRID)
  1321)       call DMDAGlobalToNaturalEnd(dm_ptr%dm,global_vec,INSERT_VALUES,natural_vec, &
  1322)                                   ierr);CHKERRQ(ierr)
  1323)     case(UNSTRUCTURED_GRID)
  1324)       call VecScatterEnd(dm_ptr%ugdm%scatter_gton,global_vec,natural_vec, &
  1325)                          INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
  1326)   end select
  1327)   
  1328) end subroutine DiscretizGlobalToNaturalEnd
  1329) 
  1330) ! ************************************************************************** !
  1331) 
  1332) subroutine DiscretizNaturalToGlobalBegin(discretization,natural_vec,global_vec,dm_index)
  1333)   ! 
  1334)   ! Begins natural to global communication with DM
  1335)   ! 
  1336)   ! Author: Glenn Hammond
  1337)   ! Date: 01/12/08
  1338)   ! 
  1339) 
  1340)   implicit none
  1341)   
  1342)   type(discretization_type) :: discretization
  1343)   Vec :: natural_vec
  1344)   Vec :: global_vec
  1345)   PetscInt :: dm_index
  1346)   PetscErrorCode :: ierr
  1347)   type(dm_ptr_type), pointer :: dm_ptr
  1348)   
  1349)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1350)   
  1351)   select case(discretization%itype)
  1352)     case(STRUCTURED_GRID)
  1353)       call DMDANaturalToGlobalBegin(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec, &
  1354)                                     ierr);CHKERRQ(ierr)
  1355)     case(UNSTRUCTURED_GRID)
  1356)   end select
  1357)   
  1358) end subroutine DiscretizNaturalToGlobalBegin
  1359) 
  1360) ! ************************************************************************** !
  1361) 
  1362) subroutine DiscretizNaturalToGlobalEnd(discretization,natural_vec,global_vec,dm_index)
  1363)   ! 
  1364)   ! Ends natural to global communication with DM
  1365)   ! 
  1366)   ! Author: Glenn Hammond
  1367)   ! Date: 01/12/08
  1368)   ! 
  1369) 
  1370)   implicit none
  1371)   
  1372)   type(discretization_type) :: discretization
  1373)   Vec :: natural_vec
  1374)   Vec :: global_vec
  1375)   PetscInt :: dm_index
  1376)   PetscErrorCode :: ierr
  1377)   type(dm_ptr_type), pointer :: dm_ptr
  1378)   
  1379)   dm_ptr => DiscretizationGetDMPtrFromIndex(discretization,dm_index)
  1380)   
  1381)   select case(discretization%itype)
  1382)     case(STRUCTURED_GRID)
  1383)       call DMDANaturalToGlobalEnd(dm_ptr%dm,natural_vec,INSERT_VALUES,global_vec, &
  1384)                                   ierr);CHKERRQ(ierr)
  1385)     case(UNSTRUCTURED_GRID)
  1386)   end select
  1387)   
  1388) end subroutine DiscretizNaturalToGlobalEnd
  1389) 
  1390) ! ************************************************************************** !
  1391) 
  1392) subroutine DiscretizationUpdateTVDGhosts(discretization,global_vec, &
  1393)                                          tvd_ghost_vec)
  1394)   ! 
  1395)   ! Updates tvd extended ghost cell values
  1396)   ! 
  1397)   ! Author: Glenn Hammond
  1398)   ! Date: 02/04/12
  1399)   ! 
  1400) 
  1401)   implicit none
  1402)   
  1403)   type(discretization_type) :: discretization
  1404)   Vec :: global_vec
  1405)   Vec :: tvd_ghost_vec
  1406)   PetscInt :: dm_index
  1407)   PetscErrorCode :: ierr
  1408)   
  1409)   call VecScatterBegin(discretization%tvd_ghost_scatter,global_vec, &
  1410)                        tvd_ghost_vec,INSERT_VALUES,SCATTER_FORWARD, &
  1411)                        ierr);CHKERRQ(ierr)
  1412)   call VecScatterEnd(discretization%tvd_ghost_scatter,global_vec, &
  1413)                        tvd_ghost_vec,INSERT_VALUES,SCATTER_FORWARD, &
  1414)                      ierr);CHKERRQ(ierr)
  1415)   
  1416) end subroutine DiscretizationUpdateTVDGhosts
  1417) 
  1418) ! ************************************************************************** !
  1419) 
  1420) subroutine DiscretAOApplicationToPetsc(discretization,int_array)
  1421)   ! 
  1422)   ! Maps application ordering to petsc
  1423)   ! 
  1424)   ! Author: Glenn Hammond
  1425)   ! Date: 10/12/12
  1426)   ! 
  1427) 
  1428)   implicit none
  1429)   
  1430) #include "petsc/finclude/petscao.h"  
  1431)   
  1432)   type(discretization_type) :: discretization
  1433)   PetscInt :: int_array(:)
  1434)   PetscErrorCode :: ierr
  1435)   
  1436)   AO :: ao
  1437)   
  1438)   select case(discretization%itype)
  1439)     case(STRUCTURED_GRID)
  1440)       call DMDAGetAO(discretization%dm_1dof,ao,ierr);CHKERRQ(ierr)
  1441)     case(UNSTRUCTURED_GRID)
  1442)       ao = discretization%grid%unstructured_grid%ao_natural_to_petsc
  1443)   end select
  1444)   call AOApplicationToPetsc(ao,size(int_array),int_array,ierr);CHKERRQ(ierr)
  1445)   
  1446) end subroutine DiscretAOApplicationToPetsc
  1447) 
  1448) ! **************************************************************************** !
  1449) 
  1450) subroutine DiscretizationInputRecord(discretization)
  1451)   ! 
  1452)   ! Prints ingested grid/discretization information
  1453)   ! 
  1454)   ! Author: Jenn Frederick
  1455)   ! Date: 03/30/2016
  1456)   ! 
  1457) 
  1458)   implicit none
  1459) 
  1460)   type(discretization_type), pointer :: discretization
  1461) 
  1462)   type(grid_type), pointer :: grid
  1463)   character(len=MAXWORDLENGTH) :: word, word1, word2
  1464)   PetscInt :: id = INPUT_RECORD_UNIT
  1465)   character(len=10) :: Format, iFormat
  1466)   
  1467)   Format = '(ES14.7)'
  1468)   iFormat = '(I10)'
  1469) 
  1470)   grid => discretization%grid
  1471) 
  1472)   write(id,'(a)') ' '
  1473)   write(id,'(a)') ' '
  1474)   write(id,'(a)') '---------------------------------------------------------&
  1475)        &-----------------------'
  1476)   write(id,'(a29)',advance='no') '---------------------------: '
  1477)   write(id,'(a)') 'GRID'
  1478)   write(id,'(a29)',advance='no') 'grid type: '
  1479)   select case(grid%itype)
  1480)     case(STRUCTURED_GRID)
  1481)       write(id,'(a)') trim(grid%ctype)
  1482)       write(id,'(a29)',advance='no') ': '
  1483)       write(id,'(a)') trim(grid%structured_grid%ctype)
  1484)       write(id,'(a29)',advance='no') 'number grid cells X: '
  1485)       write(word,iFormat) grid%structured_grid%nx
  1486)       write(id,'(a)') adjustl(trim(word)) 
  1487)       write(id,'(a29)',advance='no') 'number grid cells Y: '
  1488)       write(word,iFormat) grid%structured_grid%ny
  1489)       write(id,'(a)') adjustl(trim(word)) 
  1490)       write(id,'(a29)',advance='no') 'number grid cells Z: '
  1491)       write(word,iFormat) grid%structured_grid%nz
  1492)       write(id,'(a)') adjustl(trim(word)) 
  1493)       write(id,'(a29)',advance='no') 'delta-X (m): '
  1494)       write(id,'(1p10e12.4)') grid%structured_grid%dx_global
  1495)       write(id,'(a29)',advance='no') 'delta-Y (m): '
  1496)       write(id,'(1p10e12.4)') grid%structured_grid%dy_global
  1497)       write(id,'(a29)',advance='no') 'delta-Z (m): '
  1498)       write(id,'(1p10e12.4)') grid%structured_grid%dz_global
  1499)     case(EXPLICIT_UNSTRUCTURED_GRID,IMPLICIT_UNSTRUCTURED_GRID, &
  1500)          POLYHEDRA_UNSTRUCTURED_GRID)
  1501)       write(id,'(a)') trim(grid%ctype)
  1502)   end select
  1503) 
  1504)   write(id,'(a29)',advance='no') 'bounds X: '
  1505)   write(word1,Format) grid%structured_grid%bounds(X_DIRECTION,LOWER)
  1506)   write(word2,Format) grid%structured_grid%bounds(X_DIRECTION,UPPER)
  1507)   write(id,'(a)') adjustl(trim(word1)) // ' ,' // adjustl(trim(word2)) // ' m'
  1508)   write(id,'(a29)',advance='no') 'bounds Y: '
  1509)   write(word1,Format) grid%structured_grid%bounds(Y_DIRECTION,LOWER)
  1510)   write(word2,Format) grid%structured_grid%bounds(Y_DIRECTION,UPPER)
  1511)   write(id,'(a)') adjustl(trim(word1)) // ' ,' // adjustl(trim(word2)) // ' m'
  1512)   write(id,'(a29)',advance='no') 'bounds Z: '
  1513)   write(word1,Format) grid%structured_grid%bounds(Z_DIRECTION,LOWER)
  1514)   write(word2,Format) grid%structured_grid%bounds(Z_DIRECTION,UPPER)
  1515)   write(id,'(a)') adjustl(trim(word1)) // ' ,' // adjustl(trim(word2)) // ' m'
  1516) 
  1517)   write(id,'(a29)',advance='no') 'global origin: '
  1518)   write(word,Format) discretization%origin_global(X_DIRECTION)
  1519)   write(id,'(a)') '(x) ' // adjustl(trim(word)) // ' m'
  1520)   write(id,'(a29)',advance='no') ': '
  1521)   write(word,Format) discretization%origin_global(Y_DIRECTION)
  1522)   write(id,'(a)') '(y) ' // adjustl(trim(word)) // ' m'
  1523)   write(id,'(a29)',advance='no') ': '
  1524)   write(word,Format) discretization%origin_global(Z_DIRECTION)
  1525)   write(id,'(a)') '(z) ' // adjustl(trim(word)) // ' m'
  1526) 
  1527) end subroutine DiscretizationInputRecord
  1528) 
  1529) ! ************************************************************************** !
  1530) 
  1531) subroutine DiscretizationPrintInfo(discretization,grid,option)
  1532)   ! 
  1533)   ! Deallocates a discretization
  1534)   ! 
  1535)   ! Author: Glenn Hammond
  1536)   ! Date: 11/01/07
  1537)   ! 
  1538)   use Option_module
  1539)   use Grid_module
  1540)   
  1541)   implicit none
  1542)   
  1543)   type(discretization_type) :: discretization
  1544)   type(option_type) :: option
  1545)   
  1546)   type(grid_type), pointer :: grid
  1547)   
  1548)   grid => discretization%grid
  1549)   
  1550)   select case(discretization%itype)
  1551)     case(STRUCTURED_GRID)
  1552)       if (OptionPrintToScreen(option)) then
  1553)         write(*,'(/," Requested processors and decomposition = ", &
  1554)                  & i5,", npx,y,z= ",3i4)') &
  1555)             option%mycommsize,grid%structured_grid%npx, &
  1556)             grid%structured_grid%npy,grid%structured_grid%npz
  1557)         write(*,'(" Actual decomposition: npx,y,z= ",3i4,/)') &
  1558)             grid%structured_grid%npx_final,grid%structured_grid%npy_final, &
  1559)             grid%structured_grid%npz_final
  1560)       endif
  1561)       if (OptionPrintToScreen(option)) then
  1562)         write(option%fid_out,'(/," Requested processors and decomposition = ", &
  1563)                              & i5,", npx,y,z= ",3i4)') &
  1564)             option%mycommsize,grid%structured_grid%npx,grid%structured_grid%npy, &
  1565)             grid%structured_grid%npz
  1566)         write(option%fid_out,'(" Actual decomposition: npx,y,z= ",3i4,/)') &
  1567)             grid%structured_grid%npx_final,grid%structured_grid%npy_final, &
  1568)             grid%structured_grid%npz_final
  1569)       endif
  1570)   end select
  1571)   
  1572) end subroutine DiscretizationPrintInfo
  1573) 
  1574) ! ************************************************************************** !
  1575) 
  1576) subroutine DiscretizationDestroy(discretization)
  1577)   ! 
  1578)   ! Deallocates a discretization
  1579)   ! 
  1580)   ! Author: Glenn Hammond
  1581)   ! Date: 11/01/07
  1582)   ! 
  1583) 
  1584)   implicit none
  1585)   
  1586)   type(discretization_type), pointer :: discretization
  1587)   
  1588)   PetscErrorCode :: ierr
  1589)   PetscInt :: i
  1590)     
  1591)   if (.not.associated(discretization)) return
  1592)       
  1593)   select case(discretization%itype)
  1594)     case(STRUCTURED_GRID)
  1595)       if (discretization%dm_1dof%dm /= 0) then
  1596)         call DMDestroy(discretization%dm_1dof%dm,ierr);CHKERRQ(ierr)
  1597)       endif
  1598)       discretization%dm_1dof%dm = 0
  1599)       if (discretization%dm_nflowdof%dm /= 0) then
  1600)         call DMDestroy(discretization%dm_nflowdof%dm,ierr);CHKERRQ(ierr)
  1601)       endif
  1602)       discretization%dm_nflowdof%dm = 0
  1603)       if (discretization%dm_ntrandof%dm /= 0) then
  1604)         call DMDestroy(discretization%dm_ntrandof%dm,ierr);CHKERRQ(ierr)
  1605)       endif
  1606)       discretization%dm_n_stress_strain_dof%dm = 0
  1607)       if (discretization%dm_nflowdof%dm /= 0) then
  1608)         call DMDestroy(discretization%dm_n_stress_strain_dof%dm, &
  1609)                        ierr);CHKERRQ(ierr)
  1610)       endif
  1611)       discretization%dm_n_stress_strain_dof%dm = 0
  1612)       if (associated(discretization%dmc_nflowdof)) then
  1613)         do i=1,size(discretization%dmc_nflowdof)
  1614)           call DMDestroy(discretization%dmc_nflowdof(i)%dm,ierr);CHKERRQ(ierr)
  1615)         enddo
  1616)         deallocate(discretization%dmc_nflowdof)
  1617)         nullify(discretization%dmc_nflowdof)
  1618)       endif
  1619)       if (associated(discretization%dmc_ntrandof)) then
  1620)         do i=1,size(discretization%dmc_ntrandof)
  1621)           call DMDestroy(discretization%dmc_ntrandof(i)%dm,ierr);CHKERRQ(ierr)
  1622)         enddo
  1623)         deallocate(discretization%dmc_ntrandof)
  1624)         nullify(discretization%dmc_ntrandof)
  1625)       endif
  1626)     case(UNSTRUCTURED_GRID)
  1627)       if (associated(discretization%dm_1dof%ugdm)) &
  1628)         call UGridDMDestroy(discretization%dm_1dof%ugdm)
  1629)       if (associated(discretization%dm_nflowdof%ugdm)) &
  1630)         call UGridDMDestroy(discretization%dm_nflowdof%ugdm)
  1631)       if (associated(discretization%dm_ntrandof%ugdm)) &
  1632)         call UGridDMDestroy(discretization%dm_ntrandof%ugdm)
  1633)       if (associated(discretization%dm_n_stress_strain_dof%ugdm)) &
  1634)         call UGridDMDestroy(discretization%dm_n_stress_strain_dof%ugdm)
  1635) 
  1636)   end select
  1637)   if (associated(discretization%dm_1dof)) &
  1638)     deallocate(discretization%dm_1dof)
  1639)   nullify(discretization%dm_1dof)
  1640)   if (associated(discretization%dm_nflowdof)) &
  1641)     deallocate(discretization%dm_nflowdof)
  1642)   nullify(discretization%dm_nflowdof)
  1643)   if (associated(discretization%dm_ntrandof)) &
  1644)     deallocate(discretization%dm_ntrandof)
  1645)   nullify(discretization%dm_ntrandof)
  1646)   if (associated(discretization%dm_n_stress_strain_dof)) &
  1647)     deallocate(discretization%dm_n_stress_strain_dof)
  1648)   nullify(discretization%dm_n_stress_strain_dof)
  1649) 
  1650) 
  1651)   if (discretization%tvd_ghost_scatter /= 0) &
  1652)     call VecScatterDestroy(discretization%tvd_ghost_scatter)
  1653)   
  1654)   call GridDestroy(discretization%grid)
  1655)   
  1656)   deallocate(discretization)
  1657)   nullify(discretization)
  1658)   
  1659) end subroutine DiscretizationDestroy
  1660)  
  1661) end module Discretization_module

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