pmc_geomechanics.F90 coverage: 100.00 %func 79.03 %block
1) module PMC_Geomechanics_class
2)
3) use PMC_Base_class
4) use Realization_Subsurface_class
5) use Geomechanics_Realization_class
6) use PFLOTRAN_Constants_module
7)
8) implicit none
9)
10) #include "petsc/finclude/petscsys.h"
11)
12) private
13)
14) type, public, extends(pmc_base_type) :: pmc_geomechanics_type
15) class(realization_subsurface_type), pointer :: subsurf_realization
16) class(realization_geomech_type), pointer :: geomech_realization
17) contains
18) procedure, public :: Init => PMCGeomechanicsInit
19) procedure, public :: RunToTime => PMCGeomechanicsRunToTime
20) procedure, public :: GetAuxData => PMCGeomechanicsGetAuxData
21) procedure, public :: SetAuxData => PMCGeomechanicsSetAuxData
22) procedure, public :: Destroy => PMCGeomechanicsDestroy
23) end type pmc_geomechanics_type
24)
25) public :: PMCGeomechanicsCreate
26)
27) contains
28)
29) ! ************************************************************************** !
30)
31) function PMCGeomechanicsCreate()
32) !
33) ! This routine allocates and initializes a new object.
34) !
35) ! Author: Gautam Bisht, LBNL
36) ! Date: 01/01/14
37) !
38)
39) implicit none
40)
41) class(pmc_geomechanics_type), pointer :: PMCGeomechanicsCreate
42)
43) class(pmc_geomechanics_type), pointer :: pmc
44)
45) #ifdef DEBUG
46) print *, 'PMCGeomechanicsCreate%Create()'
47) #endif
48)
49) allocate(pmc)
50) call pmc%Init()
51)
52) PMCGeomechanicsCreate => pmc
53)
54) end function PMCGeomechanicsCreate
55)
56) ! ************************************************************************** !
57)
58) subroutine PMCGeomechanicsInit(this)
59) !
60) ! This routine initializes a new process model coupler object.
61) !
62) ! Author: Gautam Bisht, LBNL
63) ! Date: 01/01/14
64) !
65)
66) implicit none
67)
68) class(pmc_geomechanics_type) :: this
69)
70) #ifdef DEBUG
71) print *, 'PMCGeomechanics%Init()'
72) #endif
73)
74) call PMCBaseInit(this)
75) nullify(this%subsurf_realization)
76) nullify(this%geomech_realization)
77)
78) end subroutine PMCGeomechanicsInit
79)
80) ! ************************************************************************** !
81)
82) recursive subroutine PMCGeomechanicsRunToTime(this,sync_time,stop_flag)
83) !
84) ! This routine runs the geomechanics simulation.
85) !
86) ! Author: Gautam Bisht, LBNL
87) ! Date: 01/01/14
88) !
89)
90) use Timestepper_Base_class
91) use Option_module
92) use PM_Base_class
93) use Output_Geomechanics_module
94)
95) implicit none
96)
97) class(pmc_geomechanics_type), target :: this
98) PetscReal :: sync_time
99) PetscInt :: stop_flag
100) PetscInt :: local_stop_flag
101) PetscBool :: snapshot_plot_flag
102) PetscBool :: observation_plot_flag
103) PetscBool :: massbal_plot_flag
104) PetscBool :: checkpoint_flag
105)
106) class(pm_base_type), pointer :: cur_pm
107)
108) this%option%io_buffer = trim(this%name) // ':' // trim(this%pm_list%name)
109) call printVerboseMsg(this%option)
110)
111) ! Get data of other process-model
112) call this%GetAuxData()
113)
114) local_stop_flag = 0
115)
116) call SetOutputFlags(this)
117) snapshot_plot_flag = PETSC_FALSE
118) observation_plot_flag = PETSC_FALSE
119) massbal_plot_flag = PETSC_FALSE
120)
121) call this%timestepper%SetTargetTime(sync_time,this%option,local_stop_flag, &
122) snapshot_plot_flag, &
123) observation_plot_flag, &
124) massbal_plot_flag,checkpoint_flag)
125) call this%timestepper%StepDT(this%pm_list,local_stop_flag)
126)
127) ! Check if it is initial solve
128) if (this%timestepper%steps == 1) then
129) this%option%geomech_initial = PETSC_TRUE
130) endif
131)
132) ! Have to loop over all process models coupled in this object and update
133) ! the time step size. Still need code to force all process models to
134) ! use the same time step size if tightly or iteratively coupled.
135) cur_pm => this%pm_list
136) do
137) if (.not.associated(cur_pm)) exit
138) ! have to update option%time for conditions
139) this%option%time = this%timestepper%target_time
140) call cur_pm%UpdateSolution()
141) ! Geomechanics PM does not have an associate time
142) !call this%timestepper%UpdateDT(cur_pm)
143) cur_pm => cur_pm%next
144) enddo
145)
146) ! Run underlying process model couplers
147) if (associated(this%child)) then
148) ! Set data needed by process-model
149) call this%SetAuxData()
150) call this%child%RunToTime(this%timestepper%target_time,local_stop_flag)
151) endif
152)
153) if (this%timestepper%time_step_cut_flag) then
154) snapshot_plot_flag = PETSC_FALSE
155) endif
156) ! however, if we are using the modulus of the output_option%imod, we may
157) ! still print
158) if (mod(this%timestepper%steps,this%pm_list% &
159) output_option%periodic_snap_output_ts_imod) == 0) then
160) snapshot_plot_flag = PETSC_TRUE
161) endif
162) if (mod(this%timestepper%steps,this%pm_list%output_option% &
163) periodic_obs_output_ts_imod) == 0) then
164) observation_plot_flag = PETSC_TRUE
165) endif
166) if (mod(this%timestepper%steps,this%pm_list%output_option% &
167) periodic_msbl_output_ts_imod) == 0) then
168) massbal_plot_flag = PETSC_TRUE
169) endif
170)
171) call OutputGeomechanics(this%geomech_realization,snapshot_plot_flag, &
172) observation_plot_flag,massbal_plot_flag)
173) ! Set data needed by process-model
174) call this%SetAuxData()
175)
176) ! Run neighboring process model couplers
177) if (associated(this%peer)) then
178) call this%peer%RunToTime(sync_time,local_stop_flag)
179) endif
180)
181) stop_flag = max(stop_flag,local_stop_flag)
182)
183) end subroutine PMCGeomechanicsRunToTime
184)
185) ! ************************************************************************** !
186)
187) subroutine PMCGeomechanicsSetAuxData(this)
188) !
189) ! This routine updates data in simulation_aux that is required by other
190) ! process models.
191) !
192) ! Author: Gautam Bisht, LBNL
193) ! Date: 01/01/14
194) !
195)
196) use Option_module
197) use Grid_module
198)
199) implicit none
200)
201) #include "petsc/finclude/petscvec.h"
202) #include "petsc/finclude/petscvec.h90"
203)
204) class(pmc_geomechanics_type) :: this
205)
206) type(grid_type), pointer :: subsurf_grid
207) type(grid_type), pointer :: grid
208) PetscInt :: local_id
209) PetscScalar, pointer :: por0_p(:)
210) PetscScalar, pointer :: por_p(:)
211) PetscScalar, pointer :: strain_p(:)
212) PetscErrorCode :: ierr
213) PetscReal :: trace_epsilon
214) PetscReal :: por_new
215)
216) ! If at initialization stage, do nothing
217) if (this%timestepper%steps == 0) return
218)
219) select type(pmc => this)
220) class is(pmc_geomechanics_type)
221) if (this%option%geomech_subsurf_coupling == GEOMECH_TWO_WAY_COUPLED) then
222)
223) grid => pmc%subsurf_realization%patch%grid
224)
225) ! Save strain dataset in sim_aux%subsurf_strain
226) call VecScatterBegin(pmc%sim_aux%geomechanics_to_subsurf, &
227) pmc%geomech_realization%geomech_field%strain, &
228) pmc%sim_aux%subsurf_strain, &
229) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
230) call VecScatterEnd(pmc%sim_aux%geomechanics_to_subsurf, &
231) pmc%geomech_realization%geomech_field%strain, &
232) pmc%sim_aux%subsurf_strain, &
233) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
234)
235) ! Save stress dataset in sim_aux%subsurf_stress
236) call VecScatterBegin(pmc%sim_aux%geomechanics_to_subsurf, &
237) pmc%geomech_realization%geomech_field%stress, &
238) pmc%sim_aux%subsurf_stress, &
239) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
240) call VecScatterEnd(pmc%sim_aux%geomechanics_to_subsurf, &
241) pmc%geomech_realization%geomech_field%stress, &
242) pmc%sim_aux%subsurf_stress, &
243) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
244)
245) ! Update porosity dataset in sim_aux%subsurf_por
246) call VecGetArrayF90(pmc%sim_aux%subsurf_por0, por0_p, &
247) ierr);CHKERRQ(ierr)
248) call VecGetArrayF90(pmc%sim_aux%subsurf_por, por_p, &
249) ierr);CHKERRQ(ierr)
250) call VecGetArrayF90(pmc%sim_aux%subsurf_strain, strain_p, &
251) ierr);CHKERRQ(ierr)
252)
253) do local_id = 1, grid%nlmax
254) trace_epsilon = strain_p((local_id - 1)*SIX_INTEGER + ONE_INTEGER) + &
255) strain_p((local_id - 1)*SIX_INTEGER + TWO_INTEGER) + &
256) strain_p((local_id - 1)*SIX_INTEGER + THREE_INTEGER)
257) por_new = por0_p(local_id)/(1.d0 + (1.d0 - por0_p(local_id))*trace_epsilon)
258) por_p(local_id) = por_new
259) enddo
260)
261) call VecRestoreArrayF90(pmc%sim_aux%subsurf_por0, por0_p, &
262) ierr);CHKERRQ(ierr)
263) call VecRestoreArrayF90(pmc%sim_aux%subsurf_strain, strain_p, &
264) ierr);CHKERRQ(ierr)
265) call VecRestoreArrayF90(pmc%sim_aux%subsurf_por, por_p, &
266) ierr);CHKERRQ(ierr)
267)
268) endif
269)
270) end select
271)
272) end subroutine PMCGeomechanicsSetAuxData
273)
274) ! ************************************************************************** !
275)
276) subroutine PMCGeomechanicsGetAuxData(this)
277) !
278) ! This routine updates data for geomechanics simulation from other process
279) ! models.
280) !
281) ! Author: Gautam Bisht, LBNL
282) ! Date: 01/01/14
283) !
284)
285) use Option_module
286) use Geomechanics_Discretization_module
287) use Geomechanics_Force_module
288)
289) implicit none
290)
291) #include "petsc/finclude/petscsys.h"
292) #include "petsc/finclude/petscvec.h"
293)
294) class(pmc_geomechanics_type) :: this
295)
296) PetscErrorCode :: ierr
297)
298) select type(pmc => this)
299) class is(pmc_geomechanics_type)
300)
301) call VecScatterBegin(pmc%sim_aux%subsurf_to_geomechanics, &
302) pmc%sim_aux%subsurf_pres, &
303) pmc%geomech_realization%geomech_field%press, &
304) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
305) call VecScatterEnd(pmc%sim_aux%subsurf_to_geomechanics, &
306) pmc%sim_aux%subsurf_pres, &
307) pmc%geomech_realization%geomech_field%press, &
308) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
309)
310) call VecScatterBegin(pmc%sim_aux%subsurf_to_geomechanics, &
311) pmc%sim_aux%subsurf_temp, &
312) pmc%geomech_realization%geomech_field%temp, &
313) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
314) call VecScatterEnd(pmc%sim_aux%subsurf_to_geomechanics, &
315) pmc%sim_aux%subsurf_temp, &
316) pmc%geomech_realization%geomech_field%temp, &
317) INSERT_VALUES,SCATTER_FORWARD,ierr);CHKERRQ(ierr)
318)
319) call GeomechDiscretizationGlobalToLocal( &
320) pmc%geomech_realization%geomech_discretization, &
321) pmc%geomech_realization%geomech_field%press, &
322) pmc%geomech_realization%geomech_field%press_loc, &
323) ONEDOF)
324)
325) call GeomechDiscretizationGlobalToLocal( &
326) pmc%geomech_realization%geomech_discretization, &
327) pmc%geomech_realization%geomech_field%temp, &
328) pmc%geomech_realization%geomech_field%temp_loc, &
329) ONEDOF)
330)
331) end select
332)
333) end subroutine PMCGeomechanicsGetAuxData
334)
335) ! ************************************************************************** !
336)
337) subroutine PMCGeomechanicsStrip(this)
338) !
339) ! Deallocates members of PMC Geomechanics.
340) !
341) ! Author: Satish Karra
342) ! Date: 06/01/16
343)
344) implicit none
345)
346) class(pmc_geomechanics_type) :: this
347)
348) call PMCBaseStrip(this)
349) ! realizations destroyed elsewhere
350) nullify(this%subsurf_realization)
351) nullify(this%geomech_realization)
352)
353) end subroutine PMCGeomechanicsStrip
354)
355) ! ************************************************************************** !
356)
357) recursive subroutine PMCGeomechanicsDestroy(this)
358) !
359) ! Author: Satish Karra
360) ! Date: 06/01/16
361) !
362) use Option_module
363)
364) implicit none
365)
366) class(pmc_geomechanics_type) :: this
367)
368) #ifdef DEBUG
369) call printMsg(this%option,'PMCGeomechanics%Destroy()')
370) #endif
371)
372) if (associated(this%child)) then
373) call this%child%Destroy()
374) ! destroy does not currently destroy; it strips
375) deallocate(this%child)
376) nullify(this%child)
377) endif
378)
379) if (associated(this%peer)) then
380) call this%peer%Destroy()
381) ! destroy does not currently destroy; it strips
382) deallocate(this%peer)
383) nullify(this%peer)
384) endif
385)
386) call PMCGeomechanicsStrip(this)
387)
388) end subroutine PMCGeomechanicsDestroy
389)
390) end module PMC_Geomechanics_class