factory_hydrogeophysics.F90       coverage:  66.67 %func     77.91 %block


     1) module Factory_Hydrogeophysics_module
     2) 
     3)   use Simulation_Hydrogeophysics_class
     4)   
     5)   use PFLOTRAN_Constants_module
     6) 
     7)   implicit none
     8) 
     9)   private
    10) 
    11) #include "petsc/finclude/petscsys.h"
    12) 
    13)   public :: HydrogeophysicsInitialize
    14) 
    15) contains
    16) 
    17) ! ************************************************************************** !
    18) subroutine HydrogeophysicsInitialize(simulation)
    19)   ! 
    20)   ! Sets up hydrogeophysics simulation
    21)   ! 
    22)   ! Author: Glenn Hammond
    23)   ! Date: 06/17/13
    24)   ! 
    25) 
    26)   use Option_module
    27)   use Wrapper_Hydrogeophysics_module
    28)   use Input_Aux_module
    29)   use Simulation_Base_class 
    30)   use PM_Base_class  
    31)   use Discretization_module
    32)   use String_module
    33)   
    34)   implicit none
    35) 
    36) #include "petsc/finclude/petscvec.h"
    37) #include "petsc/finclude/petscvec.h90"
    38) #include "petsc/finclude/petscis.h"
    39) #include "petsc/finclude/petscviewer.h"
    40)   
    41)   class(simulation_hydrogeophysics_type) :: simulation
    42) 
    43)   type(option_type), pointer :: option
    44)   PetscMPIInt :: mycolor_mpi, mykey_mpi
    45)   PetscInt :: i, num_pflotran_processes, offset, num_slaves
    46)   PetscInt :: local_size
    47)   PetscBool :: option_found
    48)   PetscMPIInt :: mpi_int, process_range(3)
    49)   character(len=MAXSTRINGLENGTH) :: string
    50)   PetscErrorCode :: ierr
    51)   PetscInt, allocatable :: int_array(:)
    52)   IS :: is_natural
    53)   IS :: is_petsc
    54)   PetscInt :: istart
    55)   PetscInt :: temp_int
    56)   PetscViewer :: viewer
    57)   Vec :: pflotran_tracer_vec_mpi, pflotran_tracer_vec_seq
    58)   Vec :: pflotran_saturation_vec_mpi, pflotran_saturation_vec_seq
    59)   Vec :: pflotran_temperature_vec_mpi, pflotran_temperature_vec_seq
    60)   VecScatter :: pflotran_scatter
    61)   
    62) #ifdef PETSC_HAVE_MPIUNI
    63)   option%io_buffer = 'HYDROGEOPHYSICs not supported with MPIUNI.'
    64)   call printErrMsg(option)
    65) #else
    66) 
    67)   option => simulation%option
    68)   ! initialize PETSc Vecs to 0
    69)   pflotran_tracer_vec_mpi = 0
    70)   pflotran_tracer_vec_seq = 0
    71)   pflotran_saturation_vec_mpi = 0
    72)   pflotran_saturation_vec_seq = 0
    73)   pflotran_temperature_vec_mpi = 0
    74)   pflotran_temperature_vec_seq = 0
    75) 
    76)   string = '-num_slaves'
    77)   num_slaves = UNINITIALIZED_INTEGER
    78)   call InputGetCommandLineInt(string,i,option_found,option)
    79)   if (option_found) num_slaves = i
    80) 
    81)   ! NOTE: PETSc must already have been initialized here!
    82)   if (option%mycommsize < 3) then
    83)     option%io_buffer = 'At least 3 processes must be allocated to ' // &
    84)       'simulation in order to run hydrogeophysics.'
    85)     call printErrMsg(option)
    86)   endif
    87) 
    88)   ! store original communicator settings to be set back in HydrogeophysicsStrip
    89)   ! in support of multirealization scenarios.
    90)   simulation%mycomm_save = option%mycomm
    91)   simulation%myrank_save = option%myrank
    92)   simulation%mycommsize_save = option%mycommsize
    93)   simulation%mygroup_save = option%mygroup
    94)   simulation%mygroup_id_save = option%mygroup_id
    95)   
    96)   if (Uninitialized(num_slaves)) then
    97)     num_pflotran_processes = simulation%mycommsize_save / 2
    98)     num_slaves = simulation%mycommsize_save - num_pflotran_processes - 1
    99)   else if (num_slaves <= 0) then
   100)     option%io_buffer = 'Number of slaves must be greater than zero. ' // &
   101)       'Currently set to ' // StringFormatInt(num_slaves) // '.'
   102)     call printErrMsg(option)
   103)   else
   104)     if (num_slaves > simulation%mycommsize_save - 2) then
   105)       option%io_buffer = 'Too many slave processes allocated to ' // &
   106)         'simulation: ' // StringFormatInt(num_slaves)
   107)       call printErrMsg(option)
   108)     endif
   109)     num_pflotran_processes = simulation%mycommsize_save - num_slaves - 1
   110)   endif
   111)   
   112)   write(option%io_buffer,*) 'Number of E4D processes: ', &
   113)     StringFormatInt(num_slaves+1)
   114)   call printMsg(option)
   115)   write(option%io_buffer,*) 'Number of PFLOTRAN processes: ', &
   116)     StringFormatInt(num_pflotran_processes)
   117)   call printMsg(option)
   118)   
   119)   ! split the communicator
   120)   option%mygroup_id = 0
   121)   offset = 0
   122)   if (simulation%myrank_save > num_pflotran_processes-1) then
   123)     option%mygroup_id = 1
   124)     offset = num_pflotran_processes
   125)   endif
   126) 
   127)   if (option%mygroup_id == 0) then
   128)     simulation%pflotran_process = PETSC_TRUE
   129)   endif
   130) 
   131)   mycolor_mpi = option%mygroup_id
   132)   mykey_mpi = simulation%myrank_save - offset
   133)   call MPI_Comm_split(simulation%mycomm_save,mycolor_mpi,mykey_mpi, &
   134)                       option%mycomm,ierr)
   135)   call MPI_Comm_group(option%mycomm,option%mygroup,ierr)  
   136)   call MPI_Comm_rank(option%mycomm,option%myrank, ierr)
   137)   call MPI_Comm_size(option%mycomm,option%mycommsize,ierr)
   138)   
   139)   ! create group shared by both master processes
   140)   call MPI_Comm_group(simulation%mycomm_save,simulation%mygroup_save,ierr)
   141)   call MPI_Group_size(option%mygroup,mpi_int,ierr)
   142) !print *, 1, simulation%myrank_save, simulation%mygroup_save, simulation%mycomm_save
   143) !print *, 2, simulation%myrank_save, option%myrank, option%mygroup, option%mycomm, mpi_int
   144)   mpi_int = 1
   145)   process_range(1) = 0
   146)   process_range(2) = num_pflotran_processes ! includes e4d master due to 
   147)   process_range(3) = 1                        ! zero-based indexing
   148)   call MPI_Group_range_incl(simulation%mygroup_save,mpi_int,process_range, &
   149)                             simulation%pf_e4d_scatter_grp,ierr)
   150) !print *, 3, simulation%myrank_save, simulation%pf_e4d_scatter_grp
   151)   call MPI_Comm_create(simulation%mycomm_save,simulation%pf_e4d_scatter_grp, &
   152)                        simulation%pf_e4d_scatter_comm,ierr)
   153) !print *, 4, simulation%myrank_save, simulation%pf_e4d_scatter_comm, simulation%pf_e4d_master_comm
   154)   if (simulation%pf_e4d_scatter_comm /= MPI_COMM_NULL) then
   155)     call MPI_Comm_rank(simulation%pf_e4d_scatter_comm, &
   156)                        simulation%pf_e4d_scatter_rank,ierr)
   157)     call MPI_Comm_size(simulation%pf_e4d_scatter_comm, &
   158)                        simulation%pf_e4d_scatter_size,ierr)
   159) !    PFE4D_COMM = simulation%pf_e4d_scatter_comm
   160) !print *, 5, simulation%myrank_save, simulation%pf_e4d_scatter_rank, simulation%pf_e4d_scatter_size
   161)     ! remove processes between pf_master and e4d_master
   162)     mpi_int = 1
   163)     process_range(1) = 1
   164)     process_range(2) = num_pflotran_processes-1
   165)     process_range(3) = 1
   166)     ! if there are no process ranks to remove, set mpi_int to zero
   167)     if (process_range(2) - process_range(1) < 0) mpi_int = 0
   168)     call MPI_Group_range_excl(simulation%pf_e4d_scatter_grp,mpi_int, &
   169)                               process_range,simulation%pf_e4d_master_grp,ierr)
   170) !print *, 6, simulation%myrank_save, simulation%pf_e4d_master_grp, ierr
   171)     call MPI_Comm_create(simulation%pf_e4d_scatter_comm, &
   172)                          simulation%pf_e4d_master_grp, &
   173)                          simulation%pf_e4d_master_comm,ierr)
   174) !print *, 7, simulation%myrank_save, simulation%pf_e4d_master_comm
   175)     if (simulation%pf_e4d_master_comm /= MPI_COMM_NULL) then
   176)       call MPI_Comm_rank(simulation%pf_e4d_master_comm, &
   177)                          simulation%pf_e4d_master_rank,ierr)
   178)       call MPI_Comm_size(simulation%pf_e4d_master_comm, &
   179)                          simulation%pf_e4d_master_size,ierr)
   180) !      PFE4D_MASTER_COMM = simulation%pf_e4d_master_comm
   181) !print *, 8, simulation%myrank_save, simulation%pf_e4d_master_rank, simulation%pf_e4d_master_size
   182)     endif
   183)   endif
   184) 
   185) !print *, 9, simulation%myrank_save, simulation%pf_e4d_scatter_comm, simulation%pf_e4d_master_comm, MPI_COMM_NULL
   186)   
   187)   if (simulation%pflotran_process) then
   188)     call HydrogeophysicsInitPostPetsc(simulation)
   189)   else
   190)     option%io_rank = -1 ! turn off I/O from E4D processes.
   191)   endif
   192) #endif
   193) 
   194) !#define DEBUG
   195) 
   196)   !   PFLOTRAN subsurface processes      E4D master process
   197)   if (simulation%pflotran_process .or. option%myrank == 0) then 
   198)     if (simulation%pflotran_process) then
   199)       simulation%hydrogeophysics_coupler%pf_to_e4d_master_comm = &
   200)         simulation%pf_e4d_master_comm
   201)     endif
   202) 
   203)     ! create mpi Vec that includes all PFLOTRAN processes and the E4D master
   204)     call VecCreate(simulation%pf_e4d_scatter_comm,pflotran_tracer_vec_mpi, &
   205)                    ierr);CHKERRQ(ierr)
   206)     if (simulation%pflotran_process) then
   207)       call VecGetLocalSize(simulation%realization%field%work,local_size, &
   208)                            ierr);CHKERRQ(ierr)
   209)     else ! E4D master process
   210)       local_size = 0
   211)     endif
   212)     call VecSetSizes(pflotran_tracer_vec_mpi,local_size,PETSC_DECIDE, &
   213)                      ierr);CHKERRQ(ierr)
   214)     call VecSetFromOptions(pflotran_tracer_vec_mpi,ierr);CHKERRQ(ierr)
   215) 
   216)     allocate(int_array(local_size))
   217)     if (simulation%pflotran_process) then
   218)       int_array = 0
   219)       do i = 1, local_size
   220)         int_array(i) =  simulation%realization%patch%grid%nG2A( &
   221)                           simulation%realization%patch%grid%nL2G(i))
   222) 
   223)       enddo
   224)       int_array = int_array - 1
   225)     endif
   226)     call ISCreateGeneral(simulation%pf_e4d_scatter_comm,local_size,int_array, &
   227) !    call ISCreateGeneral(PETSC_COMM_SELF,local_size,int_array, &
   228)                          PETSC_COPY_VALUES,is_natural,ierr);CHKERRQ(ierr)
   229)     deallocate(int_array)
   230)     if (simulation%pflotran_process) then
   231)       call VecGetOwnershipRange(simulation%realization%field%work,istart, &
   232)                                 PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr)
   233)     endif
   234) !    istart = 0
   235)     call ISCreateStride(simulation%pf_e4d_scatter_comm,local_size,istart,1, &
   236)                         is_petsc,ierr);CHKERRQ(ierr)
   237) 
   238) #ifdef DEBUG
   239)   string = 'is_petsc.txt'
   240)   call PetscViewerASCIIOpen(simulation%pf_e4d_scatter_comm,trim(string),viewer, &
   241)                             ierr);CHKERRQ(ierr)
   242)   call ISView(is_petsc,viewer,ierr);CHKERRQ(ierr)
   243)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   244)   string = 'is_natural.txt'
   245)   call PetscViewerASCIIOpen(simulation%pf_e4d_scatter_comm,trim(string),viewer, &
   246)                             ierr);CHKERRQ(ierr)
   247) !  call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer,ierr)
   248)   call ISView(is_natural,viewer,ierr);CHKERRQ(ierr)
   249)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   250) #endif
   251) 
   252) 
   253)     ! create seq Vec on each (all PFLOTRAN processes and the E4D master)
   254)     call VecGetSize(pflotran_tracer_vec_mpi,local_size,ierr);CHKERRQ(ierr)
   255)     if (simulation%pflotran_process) local_size = 0
   256)     ! make global size the local size of pflotran_tracer_vec_seq on E4D master
   257) !    call VecCreate(PETSC_COMM_SELF,pflotran_tracer_vec_seq,ierr)
   258)     call VecCreate(simulation%pf_e4d_scatter_comm,pflotran_tracer_vec_seq, &
   259)                    ierr);CHKERRQ(ierr)
   260)     call VecSetSizes(pflotran_tracer_vec_seq,local_size,PETSC_DECIDE, &
   261)                      ierr);CHKERRQ(ierr)
   262)     call VecSetFromOptions(pflotran_tracer_vec_seq,ierr);CHKERRQ(ierr)
   263) 
   264) #ifdef DEBUG
   265)     string = 'pflotran_tracer_vec_mpi.txt'
   266)     call PetscViewerASCIIOpen(simulation%pf_e4d_scatter_comm,trim(string),viewer, &
   267)                               ierr);CHKERRQ(ierr)
   268)     call VecView(pflotran_tracer_vec_mpi,viewer,ierr);CHKERRQ(ierr)
   269)     call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   270)     string = 'pflotran_tracer_vec_seq' // trim(StringFormatInt(option%global_rank)) // '.txt'
   271)   !  call PetscViewerASCIIOpen(PETSC_COMM_SELF,trim(string),viewer,ierr)
   272)     call PetscViewerASCIIOpen(simulation%pf_e4d_scatter_comm,trim(string),viewer, &
   273)                               ierr);CHKERRQ(ierr)
   274)     call VecView(pflotran_tracer_vec_seq,viewer,ierr);CHKERRQ(ierr)
   275)     call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   276) #endif
   277)     ! create scatter between mpi Vec and local seq Vecs (only E4D master is 
   278)     ! relevant)
   279)     call VecScatterCreate(pflotran_tracer_vec_mpi,is_petsc, &
   280)                           pflotran_tracer_vec_seq,is_natural, &
   281)                           pflotran_scatter,ierr);CHKERRQ(ierr)
   282) #ifdef DEBUG
   283)     string = 'scatter.txt'
   284)     call PetscViewerASCIIOpen(simulation%pf_e4d_scatter_comm,trim(string),viewer, &
   285)                               ierr);CHKERRQ(ierr)
   286)     call VecScatterView(pflotran_scatter,viewer,ierr);CHKERRQ(ierr)
   287)     call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   288) #endif
   289)     call ISDestroy(is_natural,ierr);CHKERRQ(ierr)
   290)     call ISDestroy(is_petsc,ierr);CHKERRQ(ierr)
   291) 
   292)     ! duplicate tracer vecs into saturation version
   293)     call VecDuplicate(pflotran_tracer_vec_mpi, &
   294)                       pflotran_saturation_vec_mpi,ierr);CHKERRQ(ierr)
   295)     call VecDuplicate(pflotran_tracer_vec_seq, &
   296)                       pflotran_saturation_vec_seq,ierr);CHKERRQ(ierr)
   297)     ! have to broadcast the flow mode to the e4d master
   298)     temp_int = option%iflowmode
   299)     if (simulation%pf_e4d_master_comm /= MPI_COMM_NULL) then
   300)       call MPI_Bcast(temp_int,ONE_INTEGER_MPI,MPI_INTEGER, &
   301)                      ZERO_INTEGER_MPI,simulation%pf_e4d_master_comm,ierr)
   302)     endif
   303)     if (temp_int /= NULL_MODE .and. temp_int /= RICHARDS_MODE) then
   304)       call VecDuplicate(pflotran_tracer_vec_mpi, &
   305)                         pflotran_temperature_vec_mpi,ierr);CHKERRQ(ierr)
   306)       call VecDuplicate(pflotran_tracer_vec_seq, &
   307)                         pflotran_temperature_vec_seq,ierr);CHKERRQ(ierr)
   308)     endif
   309)   endif
   310) !print *, 'End  -----------'
   311) 
   312)   if (simulation%pflotran_process) then
   313)     ! tracer
   314)     simulation%hydrogeophysics_coupler%tracer_mpi = pflotran_tracer_vec_mpi
   315)     simulation%hydrogeophysics_coupler%tracer_seq = pflotran_tracer_vec_seq
   316)     simulation%tracer_mpi = pflotran_tracer_vec_mpi
   317)     ! saturation
   318)     simulation%hydrogeophysics_coupler%saturation_mpi = &
   319)       pflotran_saturation_vec_mpi
   320)     simulation%hydrogeophysics_coupler%saturation_seq = &
   321)       pflotran_saturation_vec_seq
   322)     simulation%saturation_mpi = pflotran_saturation_vec_mpi
   323)     ! temperature
   324)     simulation%hydrogeophysics_coupler%temperature_mpi = &
   325)       pflotran_temperature_vec_mpi
   326)     simulation%hydrogeophysics_coupler%temperature_seq = &
   327)       pflotran_temperature_vec_seq
   328)     simulation%temperature_mpi = pflotran_temperature_vec_mpi
   329)     simulation%hydrogeophysics_coupler%pf_to_e4d_scatter = pflotran_scatter
   330)     simulation%hydrogeophysics_coupler%pf_to_e4d_master_comm = &
   331)       simulation%pf_e4d_master_comm
   332)   else
   333)     call HydrogeophysicsWrapperInit(option, &
   334)                                     pflotran_tracer_vec_mpi, &
   335)                                     pflotran_tracer_vec_seq, &
   336)                                     pflotran_saturation_vec_mpi, &
   337)                                     pflotran_saturation_vec_seq, &
   338)                                     pflotran_temperature_vec_mpi, &
   339)                                     pflotran_temperature_vec_seq, &
   340)                                     pflotran_scatter, &
   341)                                     simulation%pf_e4d_master_comm)
   342)   endif
   343) #undef DEBUG
   344) 
   345) end subroutine HydrogeophysicsInitialize
   346) 
   347) ! ************************************************************************** !
   348) 
   349) subroutine HydrogeophysicsInitPostPetsc(simulation)
   350)   ! 
   351)   ! HydrogeophysicsInitializePostPetsc: Sets up hydrogeophysics simulation
   352)   ! framework after to PETSc initialization
   353)   ! 
   354)   ! Author: Glenn Hammond
   355)   ! Date: 06/17/13
   356)   ! 
   357) 
   358)   use Factory_Subsurface_module
   359)   use PMC_Hydrogeophysics_class
   360)   use Option_module
   361)   use Logging_module
   362)   
   363)   implicit none
   364)   
   365) #include "petsc/finclude/petscvec.h"
   366) #include "petsc/finclude/petscvec.h90"
   367)   
   368)   class(simulation_hydrogeophysics_type) :: simulation
   369)   
   370)   class(pmc_hydrogeophysics_type), pointer :: hydrogeophysics_coupler
   371)   character(len=MAXSTRINGLENGTH) :: string
   372)   PetscErrorCode :: ierr
   373)   
   374)   ! Init() is called in SubsurfaceInitializePostPetsc
   375)   call SubsurfaceInitializePostPetsc(simulation)
   376)   
   377)   ! add hydrogeophysics coupler to list
   378)   hydrogeophysics_coupler => PMCHydrogeophysicsCreate()
   379)   hydrogeophysics_coupler%option => simulation%option
   380)   hydrogeophysics_coupler%realization => simulation%realization
   381)   simulation%hydrogeophysics_coupler => hydrogeophysics_coupler
   382)   ! set up logging stage
   383)   string = 'Hydrogeophysics'
   384)   call LoggingCreateStage(string,hydrogeophysics_coupler%stage)
   385)   if (associated(simulation%process_model_coupler_list%child)) then
   386)     simulation%process_model_coupler_list%child%child => hydrogeophysics_coupler
   387)   else
   388)     simulation%process_model_coupler_list%child => hydrogeophysics_coupler
   389)   endif
   390) 
   391) end subroutine HydrogeophysicsInitPostPetsc
   392) 
   393) ! ************************************************************************** !
   394) 
   395) subroutine HydrogeoInitCommandLineSettings(option)
   396)   ! 
   397)   ! Initializes hydrogeophysics settings
   398)   ! 
   399)   ! Author: Glenn Hammond
   400)   ! Date: 06/17/13
   401)   ! 
   402) 
   403)   use Option_module
   404)   use Input_Aux_module
   405)   
   406)   implicit none
   407)   
   408)   type(option_type) :: option
   409)   
   410)   character(len=MAXSTRINGLENGTH) :: string
   411)   PetscBool :: option_found
   412)   PetscBool :: bool_flag
   413)   
   414) !  string = '-dummy'
   415) !  call InputGetCommandLineTruth(string,bool_flag,option_found,option)
   416) !  if (option_found) then
   417) !    option%subsurface_simulation_type = dummy
   418) !  endif
   419)   
   420) end subroutine HydrogeoInitCommandLineSettings
   421) 
   422) end module Factory_Hydrogeophysics_module

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