factory_surfsubsurface.F90       coverage:  87.50 %func     67.88 %block


     1) module Factory_Surf_Subsurf_module
     2) 
     3)   use Simulation_Surf_Subsurf_class
     4) 
     5)   use PFLOTRAN_Constants_module
     6) 
     7)   implicit none
     8) 
     9)   private
    10) 
    11) #include "petsc/finclude/petscsys.h"
    12) 
    13)   public :: SurfSubsurfaceInitialize, &
    14)             SurfSubsurfaceReadFlowPM
    15) 
    16) contains
    17) 
    18) ! ************************************************************************** !
    19) 
    20) subroutine SurfSubsurfaceInitialize(simulation)
    21)   ! 
    22)   ! This routine
    23)   ! 
    24)   ! Author: Gautam Bisht, LBNL
    25)   ! Date: 06/28/13
    26)   ! 
    27)   
    28)   implicit none
    29)   
    30)   class(simulation_surfsubsurface_type) :: simulation
    31) 
    32)   ! NOTE: PETSc must already have been initialized here!
    33)   call SurfSubsurfaceInitializePostPETSc(simulation)
    34)   
    35) end subroutine SurfSubsurfaceInitialize
    36) 
    37) ! ************************************************************************** !
    38) 
    39) subroutine SurfSubsurfaceInitializePostPETSc(simulation)
    40)   ! 
    41)   ! This routine
    42)   ! 
    43)   ! Author: Gautam Bisht, LBNL
    44)   ! Date: 06/28/13
    45)   ! 
    46) 
    47)   use Simulation_Surface_class
    48)   use Simulation_Subsurface_class
    49)   use Factory_Surface_module
    50)   use Factory_Subsurface_module
    51)   use Option_module
    52)   use Init_Common_module
    53)   use Init_Surface_module
    54)   use Surface_Flow_module
    55)   use Surface_TH_module
    56)   use Simulation_Aux_module
    57)   use PMC_Base_class
    58)   use PMC_Surface_class
    59)   use PFLOTRAN_Constants_module
    60)   use PM_Base_class
    61)   use PM_Base_Pointer_module
    62)   use PM_Surface_class
    63)   use PM_Surface_Flow_class
    64)   use PM_Surface_TH_class
    65)   use Input_Aux_module
    66)   use Realization_Subsurface_class
    67)   use String_module
    68)   use Waypoint_module
    69)   use Realization_Surface_class
    70)   use Timestepper_Surface_class
    71)   use Logging_module
    72)   use Output_Aux_module
    73)   
    74)   implicit none
    75) #include "petsc/finclude/petscvec.h"
    76) #include "petsc/finclude/petscvec.h90"
    77) 
    78)   class(simulation_surfsubsurface_type) :: simulation
    79)   
    80)   type(option_type), pointer :: option
    81)   class(realization_subsurface_type), pointer :: subsurf_realization
    82)   class(realization_surface_type), pointer :: surf_realization
    83)   class(pmc_base_type), pointer :: cur_process_model_coupler
    84)   class(pm_surface_flow_type), pointer :: pm_surface_flow
    85)   class(pm_surface_th_type), pointer :: pm_surface_th
    86)   class(pm_base_type), pointer :: cur_pm, prev_pm
    87)   class(pmc_surface_type), pointer :: pmc_surface
    88)   class(timestepper_surface_type), pointer :: timestepper
    89)   type(input_type), pointer :: input
    90)   character(len=MAXSTRINGLENGTH) :: string
    91)   PetscInt :: init_status
    92)   VecScatter :: vscat_surf_to_subsurf
    93)   VecScatter :: vscat_subsurf_to_surf
    94)   Vec :: vec_subsurf_pres
    95)   Vec :: vec_subsurf_pres_top_bc
    96)   Vec :: vec_surf_head
    97)   type(waypoint_type), pointer :: waypoint
    98)   PetscErrorCode :: ierr
    99)   
   100)   option => simulation%option
   101)   ! process command line arguments specific to subsurface
   102)   !call SurfSubsurfInitCommandLineSettings(option)
   103)   
   104)   ! we need to remove the surface pm from the list while leaving the
   105)   ! the subsurface pm linkage intact
   106)   nullify(prev_pm)
   107)   nullify(pm_surface_flow)
   108)   nullify(pm_surface_th)
   109) 
   110)   cur_pm => simulation%process_model_list
   111)   do
   112)     if (.not.associated(cur_pm)) exit
   113)     select type(cur_pm)
   114)       class is(pm_surface_th_type)
   115)         pm_surface_th => cur_pm
   116)         if (associated(prev_pm)) then
   117)           prev_pm%next => cur_pm%next
   118)         else
   119)           simulation%process_model_list => cur_pm%next
   120)         endif
   121)       class is(pm_surface_flow_type)
   122)         pm_surface_flow => cur_pm
   123)         if (associated(prev_pm)) then
   124)           prev_pm%next => cur_pm%next
   125)         else
   126)           simulation%process_model_list => cur_pm%next
   127)         endif
   128)       class default
   129)     end select
   130)     prev_pm => cur_pm
   131)     cur_pm => cur_pm%next
   132)   enddo
   133)   call SubsurfaceInitializePostPetsc(simulation)
   134)   ! in SubsurfaceInitializePostPetsc, the first pmc in the list is set as
   135)   ! the master, we need to negate this setting
   136)   simulation%process_model_coupler_list%is_master = PETSC_FALSE
   137) 
   138)   if (associated(pm_surface_flow) .or. associated(pm_surface_th)) then
   139)     simulation%surf_realization => RealizSurfCreate(option)
   140)     surf_realization => simulation%surf_realization
   141)     surf_realization%output_option => OutputOptionDuplicate(simulation%output_option)
   142)     nullify(surf_realization%output_option%output_snap_variable_list)
   143)     nullify(surf_realization%output_option%output_obs_variable_list)
   144)     surf_realization%output_option%output_snap_variable_list => OutputVariableListCreate()
   145)     surf_realization%output_option%output_obs_variable_list => OutputVariableListCreate()
   146)     subsurf_realization => simulation%realization
   147)     surf_realization%input => InputCreate(IN_UNIT,option%input_filename,option)
   148)     surf_realization%subsurf_filename = &
   149)       subsurf_realization%discretization%filename
   150)     call SurfaceInitReadRequiredCards(simulation%surf_realization)
   151)   
   152)     call setSurfaceFlowMode(option)
   153)   
   154)     if (associated(pm_surface_flow)) then
   155)       pm_surface_flow%output_option => simulation%output_option
   156)       pmc_surface => PMCSurfaceCreate()
   157)       pmc_surface%name = 'PMCSurface'
   158)       simulation%surf_flow_process_model_coupler => pmc_surface
   159)       pmc_surface%option => option
   160)       pmc_surface%checkpoint_option => simulation%checkpoint_option
   161)       pmc_surface%waypoint_list => simulation%waypoint_list_surfsubsurface
   162)       pmc_surface%pm_list => pm_surface_flow
   163)       pmc_surface%pm_ptr%pm => pm_surface_flow
   164)       pmc_surface%surf_realization => simulation%surf_realization
   165)       pmc_surface%subsurf_realization => simulation%realization
   166)       timestepper => TimestepperSurfaceCreate()
   167)       pmc_surface%timestepper => timestepper
   168)       ! set up logging stage
   169)       string = trim(pm_surface_flow%name) // 'Surface'
   170)       call LoggingCreateStage(string,pmc_surface%stage)
   171)     endif
   172) 
   173)     if (associated(pm_surface_th)) then
   174)       pm_surface_th%output_option => simulation%output_option
   175)       pmc_surface => PMCSurfaceCreate()
   176)       pmc_surface%name = 'PMCSurface'
   177)       simulation%surf_flow_process_model_coupler => pmc_surface
   178)       pmc_surface%option => option
   179)       pmc_surface%checkpoint_option => simulation%checkpoint_option
   180)       pmc_surface%waypoint_list => simulation%waypoint_list_surfsubsurface
   181)       pmc_surface%pm_list => pm_surface_th
   182)       pmc_surface%pm_ptr%pm => pm_surface_th
   183)       pmc_surface%surf_realization => simulation%surf_realization
   184)       pmc_surface%subsurf_realization => simulation%realization
   185)       timestepper => TimestepperSurfaceCreate()
   186)       pmc_surface%timestepper => timestepper
   187)       ! set up logging stage
   188)       string = trim(pm_surface_th%name) // 'Surface'
   189)       call LoggingCreateStage(string,pmc_surface%stage)
   190)     endif
   191)     
   192)     input => InputCreate(IN_UNIT,option%input_filename,option)    
   193)     string = 'SURFACE_FLOW'
   194)     call InputFindStringInFile(input,option,string)
   195)     call InputFindStringErrorMsg(input,option,string)  
   196)     call SurfaceReadInput(surf_realization,timestepper%solver, &
   197)                           simulation%waypoint_list_surfsubsurface,input)
   198)     ! Add first waypoint
   199)     waypoint => WaypointCreate()
   200)     waypoint%time = 0.d0
   201)     call WaypointInsertInList(waypoint,simulation%waypoint_list_surfsubsurface)
   202)   
   203)     ! Add final_time waypoint to surface_realization
   204)     waypoint => WaypointCreate()
   205)     waypoint%final = PETSC_TRUE
   206)     waypoint%time = simulation%waypoint_list_subsurface%last%time
   207)     waypoint%print_snap_output = PETSC_TRUE
   208)     call WaypointInsertInList(waypoint,simulation%waypoint_list_surfsubsurface)   
   209)     ! merge in outer waypoints (e.g. checkpoint times)
   210)     call WaypointListCopyAndMerge(simulation%waypoint_list_surfsubsurface, &
   211)                                   simulation%waypoint_list_subsurface,option)
   212)     call WaypointListCopyAndMerge(simulation%waypoint_list_surfsubsurface, &
   213)                                   simulation%waypoint_list_outer,option)
   214)     call InitSurfaceSetupRealization(surf_realization,subsurf_realization, &
   215)                                      simulation%waypoint_list_surfsubsurface)
   216)     call InitCommonAddOutputWaypoints(option,simulation%output_option, &
   217)                                       simulation%waypoint_list_surfsubsurface)
   218)     ! fill in holes in waypoint data
   219)     call WaypointListFillIn(simulation%waypoint_list_surfsubsurface,option)
   220)     call WaypointListRemoveExtraWaypnts(simulation% &
   221)                                             waypoint_list_surfsubsurface,option)
   222)       
   223)     if (associated(simulation%surf_flow_process_model_coupler)) then
   224)       if (associated(simulation%surf_flow_process_model_coupler% &
   225)                      timestepper)) then
   226)         simulation%surf_flow_process_model_coupler%timestepper%cur_waypoint => &
   227)           simulation%waypoint_list_surfsubsurface%first
   228)       endif
   229)     endif
   230)     call InitSurfaceSetupSolvers(surf_realization,timestepper%solver, &
   231)                              simulation%waypoint_list_surfsubsurface%last%time)
   232) 
   233)     pmc_surface%timestepper%dt_init = surf_realization%dt_init
   234)     pmc_surface%timestepper%dt_max = surf_realization%dt_max
   235)     option%surf_subsurf_coupling_flow_dt = surf_realization%dt_coupling
   236)     option%surf_flow_dt=pmc_surface%timestepper%dt_init
   237) 
   238)     if (associated(pm_surface_flow)) then
   239)       call pm_surface_flow%PMSurfaceSetRealization(surf_realization)
   240)       call pm_surface_flow%Setup()
   241)       call TSSetRHSFunction(timestepper%solver%ts, &
   242)                             pm_surface_flow%residual_vec, &
   243)                             PMRHSFunction, &
   244)                             pmc_surface%pm_ptr, &
   245)                             ierr);CHKERRQ(ierr)
   246)     endif
   247) 
   248)     if (associated(pm_surface_th)) then
   249)       call pm_surface_th%PMSurfaceSetRealization(surf_realization)
   250)       call pm_surface_th%Setup()
   251)       call TSSetRHSFunction(timestepper%solver%ts, &
   252)                             pm_surface_th%residual_vec, &
   253)                             PMRHSFunction, &
   254)                             pmc_surface%pm_ptr, &
   255)                             ierr);CHKERRQ(ierr)
   256)     endif
   257) 
   258)     timestepper%dt = option%surf_flow_dt
   259)     
   260)     nullify(simulation%process_model_coupler_list)
   261)   endif
   262) 
   263)   ! sim_aux: Create PETSc Vectors and VectorScatters
   264)   if (option%surf_flow_on .and. option%subsurf_surf_coupling /= DECOUPLED) then
   265) 
   266)     call SurfSubsurfCreateSubsurfVecs(subsurf_realization, option, &
   267)                                       vec_subsurf_pres, vec_subsurf_pres_top_bc)
   268)     call SimAuxCopySubsurfVec(simulation%sim_aux, vec_subsurf_pres)
   269)     call SimAuxCopySubsurfTopBCVec(simulation%sim_aux, vec_subsurf_pres_top_bc)
   270)     call VecDestroy(vec_subsurf_pres, ierr);CHKERRQ(ierr)
   271)     call VecDestroy(vec_subsurf_pres_top_bc, ierr);CHKERRQ(ierr)
   272) 
   273)     call SurfSubsurfCreateSurfVecs(surf_realization, option, &
   274)                                    vec_surf_head)
   275)     call SimAuxCopySurfVec(simulation%sim_aux, vec_surf_head)
   276)     call VecDestroy(vec_surf_head, ierr);CHKERRQ(ierr)
   277) 
   278)     call SurfSubsurfCreateSurfSubSurfVScats(subsurf_realization, &
   279)           surf_realization, vscat_surf_to_subsurf, vscat_subsurf_to_surf)
   280)     call SimAuxCopyVecScatter(simulation%sim_aux, vscat_surf_to_subsurf, SURF_TO_SUBSURF)
   281)     call SimAuxCopyVecScatter(simulation%sim_aux, vscat_subsurf_to_surf, SUBSURF_TO_SURF)
   282)     call VecScatterDestroy(vscat_surf_to_subsurf, ierr);CHKERRQ(ierr)
   283)     call VecScatterDestroy(vscat_surf_to_subsurf, ierr);CHKERRQ(ierr)
   284)   endif
   285) 
   286)   ! sim_aux: Set pointer
   287)   simulation%flow_process_model_coupler%sim_aux => simulation%sim_aux
   288)   if (associated(simulation%rt_process_model_coupler)) &
   289)     simulation%rt_process_model_coupler%sim_aux => simulation%sim_aux
   290)   if (option%surf_flow_on .and. &
   291)      associated(simulation%surf_flow_process_model_coupler)) &
   292)     simulation%surf_flow_process_model_coupler%sim_aux => simulation%sim_aux
   293) 
   294)   ! set surface flow as master
   295)   simulation%surf_flow_process_model_coupler%is_master = PETSC_TRUE
   296)   ! link surface flow and master
   297)   simulation%process_model_coupler_list => &
   298)     simulation%surf_flow_process_model_coupler
   299)   ! link subsurface flow as peer
   300)   simulation%process_model_coupler_list%peer => &
   301)     simulation%flow_process_model_coupler
   302) 
   303)   ! Set data in sim_aux
   304)   cur_process_model_coupler => simulation%process_model_coupler_list
   305)   call cur_process_model_coupler%SetAuxData()
   306)   if (associated(cur_process_model_coupler%peer)) then
   307)     cur_process_model_coupler => cur_process_model_coupler%peer
   308)     call cur_process_model_coupler%SetAuxData()
   309)   endif
   310) 
   311) end subroutine SurfSubsurfaceInitializePostPETSc
   312) 
   313) ! ************************************************************************** !
   314) 
   315) subroutine SurfSubsurfaceReadFlowPM(input, option, pm)
   316)   ! 
   317)   ! Author: Gautam Bisht, LBNL
   318)   ! Date: 11/29/15
   319)   !
   320)   use Input_Aux_module
   321)   use Option_module
   322)   use String_module
   323) 
   324)   use PMC_Base_class
   325)   use PM_Base_class
   326)   use PM_Surface_Flow_class
   327)   use PM_Surface_TH_class
   328)   use Init_Common_module
   329) 
   330)   implicit none
   331) 
   332)   type(input_type), pointer :: input
   333)   type(option_type), pointer :: option
   334)   class(pm_base_type), pointer :: pm
   335) 
   336)   character(len=MAXWORDLENGTH) :: word
   337)   character(len=MAXSTRINGLENGTH) :: error_string
   338) 
   339)   error_string = 'SIMULATION,PROCESS_MODELS,SURFACE_SUBSURFACE'
   340) 
   341)   option%surf_flow_on = PETSC_TRUE
   342) 
   343)   nullify(pm)
   344)   word = ''
   345)   do
   346)     call InputReadPflotranString(input,option)
   347)     if (InputCheckExit(input,option)) exit
   348)     call InputReadWord(input,option,word,PETSC_FALSE)
   349)     call StringToUpper(word)
   350)     select case(word)
   351)       case('MODE')
   352)         call InputReadWord(input,option,word,PETSC_FALSE)
   353)         call InputErrorMsg(input,option,'mode',error_string)
   354)         call StringToUpper(word)
   355)         select case(word)
   356)           case('RICHARDS')
   357)             pm => PMSurfaceFlowCreate()
   358)           case('TH')
   359)             pm => PMSurfaceTHCreate()
   360)           case default
   361)             error_string = trim(error_string) // ',MODE'
   362)             call InputKeywordUnrecognized(word,error_string,option)
   363)         end select
   364)         pm%option => option
   365)         exit
   366)       case default
   367)         error_string = trim(error_string) // ',SURFACE_FLOW'
   368)         call InputKeywordUnrecognized(word,error_string,option)
   369)     end select
   370)   enddo
   371) 
   372)   if (.not.associated(pm)) then
   373)     option%io_buffer = 'A flow MODE (card) must be included in the ' // &
   374)       'SURFACE_SUBSURFACE block in ' // trim(error_string) // '.'
   375)     call printErrMsg(option)
   376)   endif
   377) 
   378) end subroutine SurfSubsurfaceReadFlowPM
   379) 
   380) ! ************************************************************************** !
   381) 
   382) subroutine SurfSubsurfCreateSurfSubSurfVScats(realization, surf_realization, &
   383)                                               surf_to_subsurf, subsurf_to_surf)
   384)   ! 
   385)   ! This routine creates VecScatter between surface-subsurface grids.
   386)   ! Algorithm:
   387)   ! - It uses a similar logic of Matrix-Vector multiplication used in
   388)   ! UGridMapSideSet() subroutine. The algorithm here is extended to use
   389)   ! Matrix-Matrix mulitplication
   390)   ! 
   391)   ! Author: Gautam Bisht,LBNL
   392)   ! Date: 08/20/13
   393)   ! 
   394) 
   395)   use Grid_module
   396)   use String_module
   397)   use Grid_Unstructured_module
   398)   use Grid_Unstructured_Aux_module
   399)   use Grid_Unstructured_Cell_module
   400)   use Realization_Subsurface_class
   401)   use Option_module
   402)   use Patch_module
   403)   use Region_module
   404)   use Realization_Surface_class
   405) 
   406)   implicit none
   407) 
   408) #include "petsc/finclude/petscvec.h"
   409) #include "petsc/finclude/petscvec.h90"
   410) #include "petsc/finclude/petscmat.h"
   411) #include "petsc/finclude/petscmat.h90"
   412) 
   413)   class(realization_subsurface_type),pointer :: realization
   414)   class(realization_surface_type),pointer :: surf_realization
   415)   VecScatter :: surf_to_subsurf
   416)   VecScatter :: subsurf_to_surf
   417) 
   418)   type(option_type),pointer :: option
   419)   type(grid_unstructured_type),pointer :: subsurf_grid
   420)   type(grid_unstructured_type),pointer :: surf_grid
   421)   type(patch_type),pointer :: cur_patch
   422)   type(region_type),pointer :: cur_region,top_region
   423)   type(region_type),pointer :: patch_region
   424) 
   425)   Mat :: Mat_vert_to_face_subsurf
   426)   Mat :: Mat_vert_to_face_subsurf_transp
   427)   Mat :: Mat_vert_to_face_surf
   428)   Mat :: Mat_vert_to_face_surf_transp
   429)   Mat :: prod
   430)   Vec :: subsurf_petsc_ids,surf_petsc_ids
   431) 
   432)   PetscViewer :: viewer
   433) 
   434)   character(len=MAXSTRINGLENGTH) :: string
   435)   PetscInt,pointer ::int_array(:)
   436)   PetscInt :: offset
   437)   PetscInt :: int_array4(4)
   438)   PetscInt :: int_array4_0(4,1)
   439)   PetscInt :: nvertices
   440)   PetscInt :: iface
   441)   PetscInt :: local_id,ii,jj
   442)   PetscInt :: cell_type
   443)   PetscInt :: ivertex,vertex_id_local
   444)   PetscReal :: real_array4(4)
   445)   PetscReal,pointer :: vec_ptr(:)
   446) 
   447)   PetscErrorCode :: ierr
   448)   PetscBool :: found
   449) 
   450)   found = PETSC_FALSE
   451) 
   452)   if (.not.associated(realization)) return
   453) 
   454)   option => realization%option
   455)   subsurf_grid => realization%discretization%grid%unstructured_grid
   456)   surf_grid    => surf_realization%discretization%grid%unstructured_grid
   457) 
   458)   ! localize the regions on each patch
   459)   cur_patch => realization%patch_list%first
   460)   do
   461)     if (.not.associated(cur_patch)) exit
   462)     cur_region => cur_patch%region_list%first
   463)       do
   464)         if (.not.associated(cur_region)) exit
   465)         if (StringCompare(cur_region%name,'top')) then
   466)           found = PETSC_TRUE
   467)           top_region => cur_region
   468)           exit
   469)         endif
   470)         cur_region => cur_region%next
   471)       enddo
   472)     cur_patch => cur_patch%next
   473)   enddo
   474) 
   475)   if (found.eqv.PETSC_FALSE) then
   476)     option%io_buffer = 'When running with -DSURFACE_FLOW need to specify ' // &
   477)       ' in the inputfile explicitly region: top '
   478)     call printErrMsg(option)
   479)   endif
   480) 
   481)   call MatCreateAIJ(option%mycomm, &
   482)                        top_region%num_cells, &
   483)                        PETSC_DECIDE, &
   484)                        PETSC_DETERMINE, &
   485)                        subsurf_grid%num_vertices_global, &
   486)                        4, &
   487)                        PETSC_NULL_INTEGER, &
   488)                        4, &
   489)                        PETSC_NULL_INTEGER, &
   490)                        Mat_vert_to_face_subsurf, &
   491)                        ierr);CHKERRQ(ierr)
   492)   call MatCreateAIJ(option%mycomm, &
   493)                        PETSC_DECIDE, &
   494)                        top_region%num_cells, &
   495)                        subsurf_grid%num_vertices_global, &
   496)                        PETSC_DETERMINE, &
   497)                        12, &
   498)                        PETSC_NULL_INTEGER, &
   499)                        12, &
   500)                        PETSC_NULL_INTEGER, &
   501)                        Mat_vert_to_face_subsurf_transp, &
   502)                        ierr);CHKERRQ(ierr)
   503)   call VecCreateMPI(option%mycomm,top_region%num_cells,PETSC_DETERMINE, &
   504)                     subsurf_petsc_ids,ierr);CHKERRQ(ierr)
   505)   call MatZeroEntries(Mat_vert_to_face_subsurf,ierr);CHKERRQ(ierr)
   506)   real_array4 = 1.d0
   507) 
   508)   call VecGetArrayF90(subsurf_petsc_ids,vec_ptr,ierr);CHKERRQ(ierr)
   509) 
   510)   offset=0
   511)   call MPI_Exscan(top_region%num_cells,offset, &
   512)                   ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
   513) 
   514)   do ii = 1,top_region%num_cells
   515)     local_id = top_region%cell_ids(ii)
   516)     vec_ptr(ii) = subsurf_grid%cell_ids_petsc(local_id)
   517)     iface    = top_region%faces(ii)
   518)     cell_type = subsurf_grid%cell_type(local_id)
   519)     !nfaces = UCellGetNFaces(cell_type,option)
   520) 
   521)     call UCellGetNFaceVertsandVerts(option,cell_type,iface,nvertices, &
   522)                                     int_array4)
   523)     ! For this matrix:
   524)     !   irow = local face id
   525)     !   icol = natural (global) vertex id
   526)     do ivertex = 1,nvertices
   527)       vertex_id_local = &
   528)         subsurf_grid%cell_vertices(int_array4(ivertex),local_id)
   529)       int_array4_0(ivertex,1) = &
   530)         subsurf_grid%vertex_ids_natural(vertex_id_local)-1
   531)     enddo
   532)     call MatSetValues(Mat_vert_to_face_subsurf, &
   533)                       1,ii-1+offset, &
   534)                       nvertices,int_array4_0, &
   535)                       real_array4, &
   536)                       INSERT_VALUES,ierr);CHKERRQ(ierr)
   537)     call MatSetValues(Mat_vert_to_face_subsurf_transp, &
   538)                       nvertices,int_array4_0, &
   539)                       1,ii-1+offset, &
   540)                       real_array4, &
   541)                       INSERT_VALUES,ierr);CHKERRQ(ierr)
   542)   enddo
   543) 
   544)   call VecRestoreArrayF90(subsurf_petsc_ids,vec_ptr,ierr);CHKERRQ(ierr)
   545) 
   546)   call MatAssemblyBegin(Mat_vert_to_face_subsurf,MAT_FINAL_ASSEMBLY, &
   547)                         ierr);CHKERRQ(ierr)
   548)   call MatAssemblyEnd(Mat_vert_to_face_subsurf,MAT_FINAL_ASSEMBLY, &
   549)                       ierr);CHKERRQ(ierr)
   550)   call MatAssemblyBegin(Mat_vert_to_face_subsurf_transp,MAT_FINAL_ASSEMBLY, &
   551)                         ierr);CHKERRQ(ierr)
   552)   call MatAssemblyEnd(Mat_vert_to_face_subsurf_transp,MAT_FINAL_ASSEMBLY, &
   553)                       ierr);CHKERRQ(ierr)
   554) 
   555) #if UGRID_DEBUG
   556)   string = 'Mat_vert_to_face_subsurf.out'
   557)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   558)   call MatView(Mat_vert_to_face_subsurf,viewer,ierr);CHKERRQ(ierr)
   559)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   560) 
   561)   string = 'Mat_vert_to_face_subsurf_transp.out'
   562)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   563)   call MatView(Mat_vert_to_face_subsurf_transp,viewer,ierr);CHKERRQ(ierr)
   564)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   565) 
   566)   string = 'subsurf_petsc_ids.out'
   567)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   568)   call VecView(subsurf_petsc_ids,viewer,ierr);CHKERRQ(ierr)
   569)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   570) #endif
   571) 
   572) 
   573)   call MatCreateAIJ(option%mycomm, &
   574)                        surf_grid%nlmax, &
   575)                        PETSC_DECIDE, &
   576)                        PETSC_DETERMINE, &
   577)                        subsurf_grid%num_vertices_global, &
   578)                        4, &
   579)                        PETSC_NULL_INTEGER, &
   580)                        4, &
   581)                        PETSC_NULL_INTEGER, &
   582)                        Mat_vert_to_face_surf, &
   583)                        ierr);CHKERRQ(ierr)
   584)   call MatCreateAIJ(option%mycomm, &
   585)                        PETSC_DECIDE, &
   586)                        surf_grid%nlmax, &
   587)                        subsurf_grid%num_vertices_global, &
   588)                        PETSC_DETERMINE, &
   589)                        12, &
   590)                        PETSC_NULL_INTEGER, &
   591)                        12, &
   592)                        PETSC_NULL_INTEGER, &
   593)                        Mat_vert_to_face_surf_transp, &
   594)                        ierr);CHKERRQ(ierr)
   595)   call VecCreateMPI(option%mycomm,surf_grid%nlmax,PETSC_DETERMINE, &
   596)                     surf_petsc_ids,ierr);CHKERRQ(ierr)
   597)   offset=0
   598)   call MPI_Exscan(surf_grid%nlmax,offset, &
   599)                   ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
   600) 
   601)   call VecGetArrayF90(surf_petsc_ids,vec_ptr,ierr);CHKERRQ(ierr)
   602) 
   603)   do local_id = 1,surf_grid%nlmax
   604)     cell_type = surf_grid%cell_type(local_id)
   605)     vec_ptr(local_id) = surf_grid%cell_ids_petsc(local_id)
   606) 
   607)     int_array4_0 = 0
   608)     nvertices = surf_grid%cell_vertices(0,local_id)
   609)     do ivertex = 1,nvertices
   610)       vertex_id_local = surf_grid%cell_vertices(ivertex,local_id)
   611)       int_array4_0(ivertex,1) = &
   612)         surf_grid%vertex_ids_natural(vertex_id_local)-1
   613)     enddo
   614)     call MatSetValues(Mat_vert_to_face_surf, &
   615)                       1,local_id-1+offset, &
   616)                       nvertices,int_array4_0, &
   617)                       real_array4, &
   618)                       INSERT_VALUES,ierr);CHKERRQ(ierr)
   619)     call MatSetValues(Mat_vert_to_face_surf_transp, &
   620)                       nvertices,int_array4_0, &
   621)                       1,local_id-1+offset, &
   622)                       real_array4, &
   623)                       INSERT_VALUES,ierr);CHKERRQ(ierr)
   624)   enddo
   625) 
   626)   call VecRestoreArrayF90(surf_petsc_ids,vec_ptr,ierr);CHKERRQ(ierr)
   627) 
   628)   call MatAssemblyBegin(Mat_vert_to_face_surf,MAT_FINAL_ASSEMBLY, &
   629)                         ierr);CHKERRQ(ierr)
   630)   call MatAssemblyEnd(Mat_vert_to_face_surf,MAT_FINAL_ASSEMBLY, &
   631)                       ierr);CHKERRQ(ierr)
   632)   call MatAssemblyBegin(Mat_vert_to_face_surf_transp,MAT_FINAL_ASSEMBLY, &
   633)                         ierr);CHKERRQ(ierr)
   634)   call MatAssemblyEnd(Mat_vert_to_face_surf_transp,MAT_FINAL_ASSEMBLY, &
   635)                       ierr);CHKERRQ(ierr)
   636) 
   637) #if UGRID_DEBUG
   638)   string = 'Mat_vert_to_face_surf.out'
   639)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   640)   call MatView(Mat_vert_to_face_surf,viewer,ierr);CHKERRQ(ierr)
   641)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   642) 
   643)   string = 'surf_petsc_ids.out'
   644)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   645)   call VecView(surf_petsc_ids,viewer,ierr);CHKERRQ(ierr)
   646)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   647) 
   648)   string = 'Mat_vert_to_face_surf_transp.out'
   649)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   650)   call MatView(Mat_vert_to_face_surf_transp,viewer,ierr);CHKERRQ(ierr)
   651)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   652) #endif
   653) 
   654)   call MatMatMult(Mat_vert_to_face_subsurf,Mat_vert_to_face_surf_transp, &
   655)                   MAT_INITIAL_MATRIX,PETSC_DEFAULT_REAL,prod, &
   656)                   ierr);CHKERRQ(ierr)
   657) 
   658) #if UGRID_DEBUG
   659)   string = 'prod.out'
   660)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   661)   call MatView(prod,viewer,ierr);CHKERRQ(ierr)
   662)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   663) #endif
   664) 
   665)   call SurfSubsurfCreateSurfSubSurfVScat(realization,surf_realization,prod, &
   666)                                      surf_petsc_ids,surf_to_subsurf)
   667)   call MatDestroy(prod,ierr);CHKERRQ(ierr)
   668) 
   669)   call MatMatMult(Mat_vert_to_face_surf,Mat_vert_to_face_subsurf_transp, &
   670)                   MAT_INITIAL_MATRIX,PETSC_DEFAULT_REAL,prod, &
   671)                   ierr);CHKERRQ(ierr)
   672) 
   673) #if UGRID_DEBUG
   674)   string = 'prod_2.out'
   675)   call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr);CHKERRQ(ierr)
   676)   call MatView(prod,viewer,ierr);CHKERRQ(ierr)
   677)   call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
   678) #endif
   679)   call SurfSubsurfCreateSurfSubSurfVScat(realization,surf_realization,prod, &
   680)                                         subsurf_petsc_ids,subsurf_to_surf)
   681) 
   682)   call MatDestroy(prod,ierr);CHKERRQ(ierr)
   683) 
   684)   call MatDestroy(Mat_vert_to_face_subsurf,ierr);CHKERRQ(ierr)
   685)   call MatDestroy(Mat_vert_to_face_subsurf_transp,ierr);CHKERRQ(ierr)
   686)   call MatDestroy(Mat_vert_to_face_surf,ierr);CHKERRQ(ierr)
   687)   call MatDestroy(Mat_vert_to_face_surf_transp,ierr);CHKERRQ(ierr)
   688) 
   689)   call VecDestroy(subsurf_petsc_ids,ierr);CHKERRQ(ierr)
   690)   call VecDestroy(surf_petsc_ids,ierr);CHKERRQ(ierr)
   691) 
   692) end subroutine SurfSubsurfCreateSurfSubSurfVScats
   693) 
   694) ! ************************************************************************** !
   695) 
   696) subroutine SurfSubsurfCreateSurfSubSurfVScat( &
   697)               realization, &       !<
   698)               surf_realization, &  !<
   699)               prod_mat, &          !< Mat-Mat-Mult matrix
   700)               source_petsc_ids, &   !< MPI-Vector containing cell ids in PETSc order
   701)               scatter &
   702)               )
   703)   ! 
   704)   ! This subroutine creates a single vector scatter context
   705)   ! 
   706)   ! Author: Gautam Bisht,LBNL
   707)   ! Date: 08/20/13
   708)   ! 
   709) 
   710)   use Grid_module
   711)   use String_module
   712)   use Grid_Unstructured_module
   713)   use Grid_Unstructured_Cell_module
   714)   use Realization_Subsurface_class
   715)   use Option_module
   716)   use Field_module
   717)   use Surface_Field_module
   718)   use Grid_Unstructured_module
   719)   use Discretization_module
   720)   use Grid_Unstructured_Aux_module
   721)   use DM_Kludge_module
   722)   use Realization_Surface_class
   723) 
   724)   implicit none
   725) 
   726) #include "petsc/finclude/petscvec.h"
   727) #include "petsc/finclude/petscvec.h90"
   728) #include "petsc/finclude/petscmat.h"
   729) #include "petsc/finclude/petscmat.h90"
   730) 
   731)   class(realization_subsurface_type),pointer :: realization
   732)   class(realization_surface_type),pointer :: surf_realization
   733)   Mat :: prod_mat
   734)   Vec :: source_petsc_ids
   735)   VecScatter :: scatter
   736) 
   737)   Mat :: prod_loc_mat
   738)   Vec :: source_loc_vec
   739)   Vec :: corr_dest_ids_vec
   740)   Vec :: corr_dest_ids_vec_ndof
   741)   Vec :: source_petsc_ids_ndof
   742)   IS :: is_tmp1,is_tmp2
   743)   IS :: is_tmp3,is_tmp4
   744)   PetscInt,pointer :: corr_v2_ids(:)
   745)   VecScatter :: scatter_ndof
   746) 
   747)   PetscViewer :: viewer
   748) 
   749)   type(option_type),pointer :: option
   750)   type(field_type),pointer :: field
   751)   type(surface_field_type),pointer :: surf_field
   752) 
   753)   type(dm_ptr_type),pointer :: dm_ptr
   754)   character(len=MAXSTRINGLENGTH) :: string
   755)   PetscInt,pointer ::int_array(:)
   756)   PetscInt :: offset
   757)   PetscInt :: int_array4(4)
   758)   PetscInt :: int_array4_0(4,1)
   759)   PetscReal :: real_array4(4)
   760)   PetscInt :: ii,jj
   761)   PetscReal,pointer :: vec_ptr(:)
   762)   PetscInt :: ivertex,cell_id,vertex_id_local
   763)   PetscReal :: max_value
   764) 
   765)   PetscInt,pointer :: ia_p(:),ja_p(:)
   766)   PetscInt :: nrow,rstart,rend,icol(1)
   767)   PetscInt :: index
   768)   PetscInt :: vertex_id
   769)   PetscOffset :: iia,jja,iicol
   770)   PetscBool :: done
   771)   PetscScalar, pointer :: aa_v(:)
   772)   PetscInt :: row, col
   773) 
   774)   PetscErrorCode :: ierr
   775)   PetscBool :: found
   776)   PetscInt :: nlocal
   777) 
   778)   option     => realization%option
   779)   field      => realization%field
   780)   surf_field => surf_realization%surf_field
   781) 
   782)   if (option%mycommsize > 1) then
   783)     ! From the MPI-Matrix get the local-matrix
   784)     call MatMPIAIJGetLocalMat(prod_mat,MAT_INITIAL_MATRIX,prod_loc_mat, &
   785)                               ierr);CHKERRQ(ierr)
   786)     ! Get i and j indices of the local-matrix
   787)     call MatGetRowIJF90(prod_loc_mat,ONE_INTEGER,PETSC_FALSE,PETSC_FALSE, &
   788)                         nrow,ia_p,ja_p,done,ierr);CHKERRQ(ierr)
   789)     ! Get values stored in the local-matrix
   790)     call MatSeqAIJGetArrayF90(prod_loc_mat,aa_v,ierr);CHKERRQ(ierr)
   791)   else
   792)     ! Get i and j indices of the local-matrix
   793)     call MatGetRowIJF90(prod_mat,ONE_INTEGER,PETSC_FALSE,PETSC_FALSE, &
   794)                         nrow,ia_p,ja_p,done,ierr);CHKERRQ(ierr)
   795)     ! Get values stored in the local-matrix
   796)     call MatSeqAIJGetArrayF90(prod_mat,aa_v,ierr);CHKERRQ(ierr)
   797)   endif
   798) 
   799)   ! For each row of the local-matrix,find the column with the largest value
   800)   allocate(corr_v2_ids(nrow))
   801)   row = 1
   802)   col = 0
   803)   do ii = 1,nrow
   804)     max_value = 0.d0
   805)     do jj = ia_p(ii),ia_p(ii + 1) - 1
   806)       if (aa_v(jj) > max_value) then
   807)         corr_v2_ids(ii) = ja_p(jj)
   808)         max_value = aa_v(jj)
   809)       endif
   810)     enddo
   811)     if (max_value<3) then
   812)       option%io_buffer = 'Atleast three vertices need to form a face'
   813)       call printErrMsg(option)
   814)     endif
   815)   enddo
   816) 
   817)   offset = 0
   818)   call MPI_Exscan(nrow,offset,ONE_INTEGER_MPI, &
   819)                   MPIU_INTEGER,MPI_SUM,option%mycomm,ierr)
   820)   allocate(int_array(nrow))
   821)   do ii = 1,nrow
   822)     int_array(ii) = (ii-1)+offset
   823)   enddo
   824)   call ISCreateGeneral(option%mycomm,nrow, &
   825)                        int_array,PETSC_COPY_VALUES,is_tmp1,ierr);CHKERRQ(ierr)
   826)   call ISCreateBlock(option%mycomm,option%nflowdof,nrow, &
   827)                      int_array,PETSC_COPY_VALUES,is_tmp3,ierr);CHKERRQ(ierr)
   828) 
   829)   do ii = 1,nrow
   830)     int_array(ii) = corr_v2_ids(ii)-1
   831)   enddo
   832)   call ISCreateGeneral(option%mycomm,nrow, &
   833)                        int_array,PETSC_COPY_VALUES,is_tmp2,ierr);CHKERRQ(ierr)
   834)   call ISCreateBlock(option%mycomm,option%nflowdof,nrow, &
   835)                      int_array,PETSC_COPY_VALUES,is_tmp4,ierr);CHKERRQ(ierr)
   836)   deallocate(int_array)
   837) 
   838)   call VecCreateMPI(option%mycomm,nrow,PETSC_DETERMINE, &
   839)                     corr_dest_ids_vec,ierr);CHKERRQ(ierr)
   840)   call VecScatterCreate(source_petsc_ids,is_tmp2,corr_dest_ids_vec,is_tmp1, &
   841)                         scatter,ierr);CHKERRQ(ierr)
   842)   call ISDestroy(is_tmp1,ierr);CHKERRQ(ierr)
   843)   call ISDestroy(is_tmp2,ierr);CHKERRQ(ierr)
   844) 
   845)   call VecScatterBegin(scatter,source_petsc_ids,corr_dest_ids_vec, &
   846)                        INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
   847)   call VecScatterEnd(scatter,source_petsc_ids,corr_dest_ids_vec, &
   848)                        INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
   849) 
   850)   call VecDestroy(corr_dest_ids_vec,ierr);CHKERRQ(ierr)
   851)   if (option%mycommsize>1) then
   852)     call MatSeqAIJRestoreArrayF90(prod_loc_mat,aa_v,ierr);CHKERRQ(ierr)
   853)     call MatDestroy(prod_loc_mat,ierr);CHKERRQ(ierr)
   854)   else
   855)     call MatSeqAIJRestoreArrayF90(prod_mat,aa_v,ierr);CHKERRQ(ierr)
   856)   endif
   857) 
   858) end subroutine SurfSubsurfCreateSurfSubSurfVScat
   859) 
   860) ! ************************************************************************** !
   861) 
   862) subroutine SurfSubsurfCreateSubsurfVecs(subsurf_realization, option, &
   863)                                         subsurf_pres, subsurf_pres_top_bc)
   864)   ! 
   865)   ! This routine creates subsurface vectors.
   866)   ! 
   867)   ! Author: Gautam Bisht,LBNL
   868)   ! Date: 08/20/13
   869)   ! 
   870) 
   871)   use Realization_Subsurface_class
   872)   use Coupler_module
   873)   use Option_module
   874)   use String_module
   875) 
   876)   implicit none
   877) 
   878)   class(realization_subsurface_type),pointer :: subsurf_realization
   879)   type(option_type),pointer :: option
   880)   Vec :: subsurf_pres
   881)   Vec :: subsurf_pres_top_bc
   882) 
   883)   type(coupler_list_type),pointer :: coupler_list
   884)   type(coupler_type),pointer :: coupler
   885) 
   886)   PetscInt :: num_conn
   887)   PetscInt :: found
   888)   PetscInt :: found_global
   889)   PetscErrorCode :: ierr
   890) 
   891)   call VecCreate(option%mycomm,subsurf_pres,ierr);CHKERRQ(ierr)
   892)   call VecSetSizes(subsurf_pres, &
   893)                    subsurf_realization%discretization%grid%nlmax, &
   894)                    PETSC_DECIDE,ierr);CHKERRQ(ierr)
   895)   call VecSetFromOptions(subsurf_pres,ierr);CHKERRQ(ierr)
   896)   call VecSet(subsurf_pres,0.d0,ierr);CHKERRQ(ierr)
   897) 
   898) #if 1
   899)   found = 0
   900)   num_conn = 0
   901)   coupler_list => subsurf_realization%patch%source_sink_list
   902)   coupler => coupler_list%first
   903)   do
   904)     if (.not.associated(coupler)) exit
   905)     if (StringCompare(coupler%name,'from_surface_ss')) then
   906)       num_conn = coupler%connection_set%num_connections
   907)       found = 1
   908)     endif
   909)     coupler => coupler%next
   910)   enddo
   911) 
   912)   call MPI_AllReduce(found,found_global,ONE_INTEGER_MPI,MPI_INTEGER,MPI_MAX, &
   913)                      option%mycomm,ierr)
   914) #endif
   915) 
   916)   !found_global = 0
   917)   if (found_global == 0) then
   918)     coupler_list => subsurf_realization%patch%boundary_condition_list
   919)     coupler => coupler_list%first
   920)     do
   921)       if (.not.associated(coupler)) exit
   922)       if (StringCompare(coupler%name,'from_surface_bc')) then
   923)         num_conn = coupler%connection_set%num_connections
   924)         found = 1
   925)       endif
   926)       coupler => coupler%next
   927)     enddo
   928)     call MPI_AllReduce(found,found_global,ONE_INTEGER_MPI,MPI_INTEGER,MPI_MAX, &
   929)                        option%mycomm,ierr)
   930)   endif
   931) 
   932)   if (found_global > 0) then
   933)     call VecCreate(option%mycomm,subsurf_pres_top_bc,ierr);CHKERRQ(ierr)
   934)     call VecSetSizes(subsurf_pres_top_bc,num_conn,PETSC_DECIDE, &
   935)                      ierr);CHKERRQ(ierr)
   936)     call VecSetFromOptions(subsurf_pres_top_bc,ierr);CHKERRQ(ierr)
   937)     call VecSet(subsurf_pres_top_bc,0.d0,ierr);CHKERRQ(ierr)
   938)   endif
   939) 
   940) end subroutine SurfSubsurfCreateSubsurfVecs
   941) 
   942) ! ************************************************************************** !
   943) 
   944) subroutine SurfSubsurfCreateSurfVecs(surf_realization,option,surf_head)
   945)   ! 
   946)   ! This routine creates surface vectors.
   947)   ! 
   948)   ! Author: Gautam Bisht,LBNL
   949)   ! Date: 08/20/13
   950)   ! 
   951) 
   952)   use Realization_Surface_class
   953)   use Option_module
   954) 
   955)   implicit none
   956) 
   957)   class(realization_surface_type),pointer :: surf_realization
   958)   type(option_type),pointer :: option
   959)   Vec :: surf_head
   960) 
   961)   PetscErrorCode :: ierr
   962) 
   963)   call VecCreate(option%mycomm,surf_head,ierr);CHKERRQ(ierr)
   964)   call VecSetSizes(surf_head,surf_realization%discretization%grid%nlmax, &
   965)                    PETSC_DECIDE,ierr);CHKERRQ(ierr)
   966)   call VecSetFromOptions(surf_head,ierr);CHKERRQ(ierr)
   967)   call VecSet(surf_head,0.d0,ierr);CHKERRQ(ierr)
   968) 
   969) end subroutine SurfSubsurfCreateSurfVecs
   970) 
   971) ! ************************************************************************** !
   972) 
   973) subroutine SurfSubsurfInitCommandLineSettings(option)
   974)   ! 
   975)   ! This routine
   976)   ! 
   977)   ! Author: Gautam Bisht, LBNL
   978)   ! Date: 07/01/13
   979)   ! 
   980) 
   981)   use Option_module
   982)   use Input_Aux_module
   983)   
   984)   implicit none
   985)   
   986)   type(option_type) :: option
   987)   
   988)   character(len=MAXSTRINGLENGTH) :: string
   989)   PetscBool :: option_found
   990)   PetscBool :: bool_flag
   991)   
   992)   string = '-multisimulation'
   993)   call InputGetCommandLineTruth(string,bool_flag,option_found,option)
   994)   if (option_found) then
   995)     option%subsurface_simulation_type = MULTISIMULATION_SIM_TYPE
   996)   endif
   997) 
   998)   string = '-stochastic'
   999)   call InputGetCommandLineTruth(string,bool_flag,option_found,option)
  1000)   if (option_found) then
  1001)     option%subsurface_simulation_type = STOCHASTIC_SIM_TYPE
  1002)   endif
  1003)   
  1004) end subroutine SurfSubsurfInitCommandLineSettings
  1005) 
  1006) end module Factory_Surf_Subsurf_module

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