surface_flow.F90 coverage: 0.00 %func 0.00 %block
1) module Surface_Flow_module
2)
3) use Global_Aux_module
4)
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) private
10)
11) #include "petsc/finclude/petscsys.h"
12)
13) #include "petsc/finclude/petscvec.h"
14) #include "petsc/finclude/petscvec.h90"
15) #include "petsc/finclude/petscmat.h"
16) #include "petsc/finclude/petscmat.h90"
17) !#include "petsc/finclude/petscsnes.h"
18) #include "petsc/finclude/petscviewer.h"
19) #include "petsc/finclude/petsclog.h"
20) #include "petsc/finclude/petscts.h"
21)
22) ! Cutoff parameters
23) PetscReal, parameter :: eps = 1.D-12
24) PetscReal, parameter :: perturbation_tolerance = 1.d-6
25)
26) public SurfaceFlowSetup, &
27) SurfaceFlowDiffusion, &
28) SurfaceFlowUpdateSolution, &
29) SurfaceFlowRHSFunction, &
30) SurfaceFlowComputeMaxDt, &
31) SurfaceFlowGetTecplotHeader, &
32) SurfaceFlowUpdateAuxVars, &
33) SurfaceFlowUpdateSurfState
34)
35) contains
36)
37) ! ************************************************************************** !
38)
39) subroutine SurfaceFlowSetup(surf_realization)
40) !
41) ! This routine sets up surfaceflow type
42) !
43) ! Author: Gautam Bisht, ORNL
44) ! Date: 05/21/12
45) !
46)
47) use Realization_Surface_class
48) use Output_Aux_module
49)
50) class(realization_surface_type) :: surf_realization
51)
52) type(output_variable_list_type), pointer :: list
53)
54) list => surf_realization%output_option%output_snap_variable_list
55) call SurfaceFlowSetPlotVariables(list)
56) list => surf_realization%output_option%output_obs_variable_list
57) call SurfaceFlowSetPlotVariables(list)
58)
59) end subroutine SurfaceFlowSetup
60)
61) ! ************************************************************************** !
62)
63) subroutine SurfaceFlowSetPlotVariables(list)
64) !
65) ! This routine adds default variables to be printed to list
66) !
67) ! Author: Gautam Bisht, LBNL
68) ! Date: 10/30/12
69) !
70)
71) use Realization_Surface_class
72) use Output_Aux_module
73) use Variables_module
74)
75) implicit none
76)
77) type(output_variable_list_type), pointer :: list
78)
79) character(len=MAXWORDLENGTH) :: name, units
80)
81) if (associated(list%first)) then
82) return
83) endif
84)
85) name = 'H'
86) units = 'm'
87) call OutputVariableAddToList(list,name,OUTPUT_GENERIC,units, &
88) SURFACE_LIQUID_HEAD)
89)
90) name = 'Material ID'
91) units = ''
92) call OutputVariableAddToList(list,name,OUTPUT_DISCRETE,units, &
93) MATERIAL_ID)
94)
95) end subroutine SurfaceFlowSetPlotVariables
96)
97) ! ************************************************************************** !
98)
99) subroutine SurfaceFlowKinematic(hw_up, &
100) mannings_up, &
101) hw_dn, &
102) mannings_dn, &
103) slope, &
104) length, &
105) option, &
106) vel, &
107) Res)
108) !
109) ! This routine computes the internal flux term for the residual under
110) ! kinematic-wave assumption.
111) !
112) ! Author: Gautam Bisht, ORNL
113) ! Date: 05/21/12
114) !
115)
116) use Option_module
117)
118) implicit none
119)
120) type(option_type) :: option
121) PetscReal :: hw_up, hw_dn
122) PetscReal :: slope
123) PetscReal :: mannings_up, mannings_dn
124) PetscReal :: length
125) PetscReal :: vel ! units: m/s
126) PetscReal :: Res(1:option%nflowdof) ! units: m^3/s
127)
128) PetscReal :: flux ! units: m^2/s
129)
130) ! initialize
131) flux = 0.d0
132) vel = 0.d0
133)
134) if (slope<0.d0) then
135) vel = sqrt(dabs(slope))/mannings_up*((hw_up)**(2.d0/3.d0))
136) flux= hw_up*vel
137) else
138) vel = -sqrt(dabs(slope))/mannings_dn*((hw_dn)**(2.d0/3.d0))
139) flux= hw_dn*vel
140) endif
141)
142) Res(1) = flux*length
143)
144) end subroutine SurfaceFlowKinematic
145)
146) ! ************************************************************************** !
147)
148) subroutine SurfaceFlowDiffusion(hw_up, &
149) zc_up, &
150) mannings_up, &
151) hw_dn, &
152) zc_dn, &
153) mannings_dn, &
154) dist, &
155) length, &
156) option, &
157) vel, &
158) Res)
159) !
160) ! This routine computes the internal flux term for the residual under
161) ! diffusion-wave assumption.
162) !
163) ! Author: Gautam Bisht, LBL
164) ! Date: 08/03/12
165) !
166)
167) use Option_module
168)
169) implicit none
170)
171) type(option_type) :: option
172) PetscReal :: hw_up, hw_dn
173) PetscReal :: zc_up, zc_dn
174) PetscReal :: head_up, head_dn
175) PetscReal :: mannings_up, mannings_dn
176) PetscReal :: dist, length
177) PetscReal :: vel ! units: m/s
178) PetscReal :: Res(1:option%nflowdof) ! units: m^3/s
179)
180) PetscReal :: flux ! units: m^2/s
181) PetscReal :: Cd
182) PetscReal :: hw_half
183) PetscReal :: mannings_half
184) PetscReal :: dhead
185)
186) ! initialize
187) flux = 0.d0
188) Cd = 1.0d0
189)
190) head_up = hw_up + zc_up
191) head_dn = hw_dn + zc_dn
192)
193) if (head_up>head_dn) then
194) mannings_half = mannings_up
195) if (hw_up>0.d0) then
196) hw_half = hw_up
197) else
198) hw_half = 0.d0
199) endif
200) else
201) mannings_half = mannings_dn
202) if (hw_dn>0.d0) then
203) hw_half = hw_dn
204) else
205) hw_half = 0.d0
206) endif
207) endif
208)
209) dhead=head_up-head_dn
210) if (abs(dhead)<eps) then
211) dhead=0.d0
212) vel = 0.d0
213) else
214) vel = (hw_half**(2.d0/3.d0))/mannings_half* &
215) dhead/(abs(dhead)**(1.d0/2.d0))* &
216) 1.d0/(dist**0.5d0)
217) endif
218)
219) flux = hw_half*vel
220) Res(1) = flux*length
221)
222) end subroutine SurfaceFlowDiffusion
223)
224) ! ************************************************************************** !
225)
226) subroutine SurfaceFlowUpdateSolution(surf_realization)
227) !
228) ! This routine updates data in module after a successful time step
229) !
230) ! Author: Gautam Bisht, ORNL
231) ! Date: 05/22/12
232) !
233)
234) use Realization_Surface_class
235) use Surface_Field_module
236)
237) implicit none
238)
239) class(realization_surface_type) :: surf_realization
240)
241) type(surface_field_type),pointer :: surf_field
242) PetscErrorCode :: ierr
243)
244) surf_field => surf_realization%surf_field
245) call VecCopy(surf_field%flow_xx,surf_field%flow_yy,ierr);CHKERRQ(ierr)
246)
247) end subroutine SurfaceFlowUpdateSolution
248)
249) ! ************************************************************************** !
250)
251) subroutine SurfaceFlowRHSFunction(ts,t,xx,ff,surf_realization,ierr)
252) !
253) ! This routine provides the function evaluation for PETSc TSSolve()
254) !
255) ! Author: Gautam Bisht, LBNL
256) ! Date: 03/07/13
257) !
258)
259) use Realization_Surface_class
260) use Surface_Field_module
261) use Patch_module
262) use Discretization_module
263) use Option_module
264) use Logging_module
265) use Connection_module
266) use Grid_module
267) use Coupler_module
268) use Surface_Field_module
269) use Debug_module
270) use Surface_Global_Aux_module
271)
272) implicit none
273)
274) TS :: ts
275) PetscReal :: t
276) Vec :: xx
277) Vec :: ff
278) class(realization_surface_type) :: surf_realization
279) PetscErrorCode :: ierr
280)
281) PetscViewer :: viewer
282)
283) type(discretization_type), pointer :: discretization
284) type(surface_field_type), pointer :: surf_field
285) type(option_type), pointer :: option
286) type(grid_type), pointer :: grid
287) type(patch_type), pointer :: patch
288) type(coupler_type), pointer :: boundary_condition
289) type(coupler_type), pointer :: source_sink
290) type(connection_set_list_type), pointer :: connection_set_list
291) type(connection_set_type), pointer :: cur_connection_set
292) type(surface_global_auxvar_type), pointer :: surf_global_auxvars(:)
293) type(surface_global_auxvar_type), pointer :: surf_global_auxvars_bc(:)
294)
295) PetscInt :: local_id_up, local_id_dn, local_id
296) PetscInt :: ghosted_id_up, ghosted_id_dn, ghosted_id
297) PetscInt :: iconn
298) PetscInt :: sum_connection
299) PetscReal :: dx, dy, dz
300) PetscReal :: dist
301) PetscReal :: vel
302) PetscReal :: slope, slope_dn
303) PetscReal :: hw_up, hw_dn ! water height [m]
304) PetscReal :: Res(surf_realization%option%nflowdof), v_darcy
305) PetscReal :: qsrc, qsrc_flow
306)
307) character(len=MAXSTRINGLENGTH) :: string,string2
308)
309) PetscReal, pointer :: ff_p(:), mannings_loc_p(:),area_p(:)
310) PetscReal, pointer :: xc(:),yc(:),zc(:)
311)
312) patch => surf_realization%patch
313) grid => patch%grid
314) option => surf_realization%option
315) surf_field => surf_realization%surf_field
316) surf_global_auxvars => patch%surf_aux%SurfaceGlobal%auxvars
317) surf_global_auxvars_bc => patch%surf_aux%SurfaceGlobal%auxvars_bc
318)
319) surf_field => surf_realization%surf_field
320) discretization => surf_realization%discretization
321) option => surf_realization%option
322) patch => surf_realization%patch
323) grid => patch%grid
324) surf_global_auxvars => patch%surf_aux%SurfaceGlobal%auxvars
325) surf_global_auxvars_bc => patch%surf_aux%SurfaceGlobal%auxvars_bc
326)
327) surf_realization%iter_count = surf_realization%iter_count+1
328) if (surf_realization%iter_count < 10) then
329) write(string2,'("00",i1)') surf_realization%iter_count
330) else if (surf_realization%iter_count < 100) then
331) write(string2,'("0",i2)') surf_realization%iter_count
332) else if (surf_realization%iter_count < 1000) then
333) write(string2,'(i3)') surf_realization%iter_count
334) else if (surf_realization%iter_count < 10000) then
335) write(string2,'(i4)') surf_realization%iter_count
336) endif
337)
338) call DiscretizationGlobalToLocal(discretization,xx,surf_field%flow_xx_loc,NFLOWDOF)
339) ! Then, update the aux vars
340) call SurfaceFlowUpdateAuxVars(surf_realization)
341)
342) call VecGetArrayF90(ff,ff_p, ierr);CHKERRQ(ierr)
343) call VecGetArrayF90(surf_field%mannings_loc,mannings_loc_p, &
344) ierr);CHKERRQ(ierr)
345) call VecGetArrayF90(surf_field%area,area_p,ierr);CHKERRQ(ierr)
346)
347) ff_p = 0.d0
348) Res = 0.d0
349)
350) xc => surf_realization%discretization%grid%x
351) yc => surf_realization%discretization%grid%y
352) zc => surf_realization%discretization%grid%z
353)
354) ! Interior Flux Terms -----------------------------------
355) connection_set_list => grid%internal_connection_set_list
356) cur_connection_set => connection_set_list%first
357) sum_connection = 0
358) do
359) if (.not.associated(cur_connection_set)) exit
360) do iconn = 1, cur_connection_set%num_connections
361) sum_connection = sum_connection + 1
362)
363) ghosted_id_up = cur_connection_set%id_up(iconn)
364) ghosted_id_dn = cur_connection_set%id_dn(iconn)
365)
366) local_id_up = grid%nG2L(ghosted_id_up)
367) local_id_dn = grid%nG2L(ghosted_id_dn)
368)
369) dx = xc(ghosted_id_dn) - xc(ghosted_id_up)
370) dy = yc(ghosted_id_dn) - yc(ghosted_id_up)
371) dz = zc(ghosted_id_dn) - zc(ghosted_id_up)
372) dist = sqrt(dx*dx + dy*dy + dz*dz)
373) slope = dz/dist
374)
375) if (surf_global_auxvars(ghosted_id_up)%head(1)<0.d0 .or. &
376) surf_global_auxvars(ghosted_id_dn)%head(1)<0.d0) then
377) write(*,*) 'In SurfaceFlowFlux: ', surf_global_auxvars(ghosted_id_up)%head(1), &
378) surf_global_auxvars(ghosted_id_dn)%head(1),ghosted_id_up,ghosted_id_dn
379) option%io_buffer='stopping: -ve head values '
380) call printErrMsg(option)
381) endif
382)
383) select case(option%surface_flow_formulation)
384) case (KINEMATIC_WAVE)
385) option%io_buffer='Explicit Surface flow not implemented for ' // &
386) 'Kinematic wave'
387) call printErrMsg(option)
388) case (DIFFUSION_WAVE)
389) call SurfaceFlowFlux(surf_global_auxvars(ghosted_id_up), &
390) zc(ghosted_id_up), &
391) mannings_loc_p(ghosted_id_up), &
392) surf_global_auxvars(ghosted_id_dn), &
393) zc(ghosted_id_dn), &
394) mannings_loc_p(ghosted_id_dn), &
395) dist, cur_connection_set%area(iconn), &
396) option,vel,Res)
397) end select
398)
399) patch%internal_velocities(1,sum_connection) = vel
400) patch%internal_flow_fluxes(1,sum_connection) = Res(1)
401)
402) if (local_id_up>0) then
403) ff_p(local_id_up) = ff_p(local_id_up) - Res(1)/area_p(local_id_up)
404) endif
405)
406) if (local_id_dn>0) then
407) ff_p(local_id_dn) = ff_p(local_id_dn) + Res(1)/area_p(local_id_dn)
408) endif
409)
410)
411) enddo
412) cur_connection_set => cur_connection_set%next
413) enddo
414)
415) ! Boundary Flux Terms -----------------------------------
416) boundary_condition => patch%boundary_condition_list%first
417) sum_connection = 0
418) do
419) if (.not.associated(boundary_condition)) exit
420)
421) cur_connection_set => boundary_condition%connection_set
422)
423) do iconn = 1, cur_connection_set%num_connections
424) sum_connection = sum_connection + 1
425)
426) local_id = cur_connection_set%id_dn(iconn)
427) ghosted_id_dn = grid%nL2G(local_id)
428)
429) dx = xc(ghosted_id_dn) - cur_connection_set%intercp(1,iconn)
430) dy = yc(ghosted_id_dn) - cur_connection_set%intercp(2,iconn)
431) dz = zc(ghosted_id_dn) - cur_connection_set%intercp(3,iconn)
432) slope_dn = dz/sqrt(dx*dx + dy*dy + dz*dz)
433)
434) call SurfaceFlowBCFlux(boundary_condition%flow_condition%itype, &
435) surf_global_auxvars_bc(sum_connection), &
436) slope_dn, &
437) mannings_loc_p(ghosted_id_dn), &
438) cur_connection_set%area(iconn), &
439) option,vel,Res)
440)
441) patch%boundary_velocities(1,sum_connection) = vel
442) patch%boundary_flow_fluxes(1,sum_connection) = Res(1)
443)
444) ff_p(local_id) = ff_p(local_id) + Res(1)/area_p(local_id)
445) enddo
446) boundary_condition => boundary_condition%next
447) enddo
448)
449) ! Source/sink terms -------------------------------------
450) source_sink => patch%source_sink_list%first
451) sum_connection = 0
452) do
453) if (.not.associated(source_sink)) exit
454)
455) qsrc_flow = 0.d0
456) if (source_sink%flow_condition%rate%itype/=HET_VOL_RATE_SS.and. &
457) source_sink%flow_condition%rate%itype/=HET_MASS_RATE_SS) &
458) qsrc_flow = source_sink%flow_condition%rate%dataset%rarray(1)
459)
460) cur_connection_set => source_sink%connection_set
461)
462) do iconn = 1, cur_connection_set%num_connections
463) sum_connection = sum_connection + 1
464) local_id = cur_connection_set%id_dn(iconn)
465) ghosted_id = grid%nL2G(local_id)
466) if (patch%imat(ghosted_id) <= 0) cycle
467)
468) select case(source_sink%flow_condition%rate%itype)
469) case(VOLUMETRIC_RATE_SS) ! assume local density for now
470) ! qsrc = m^3/sec
471) qsrc = qsrc_flow*area_p(local_id)
472) case(HET_VOL_RATE_SS)
473) ! qsrc = m^3/sec
474) qsrc = source_sink%flow_aux_real_var(ONE_INTEGER,iconn)*area_p(local_id)
475) case default
476) option%io_buffer = 'Source/Sink flow condition type not recognized'
477) call printErrMsg(option)
478) end select
479)
480) ff_p(local_id) = ff_p(local_id) + qsrc/area_p(local_id)
481) enddo
482) source_sink => source_sink%next
483) enddo
484)
485) call VecRestoreArrayF90(ff,ff_p, ierr);CHKERRQ(ierr)
486) call VecRestoreArrayF90(surf_field%mannings_loc,mannings_loc_p, &
487) ierr);CHKERRQ(ierr)
488) call VecRestoreArrayF90(surf_field%area,area_p,ierr);CHKERRQ(ierr)
489)
490) if (surf_realization%debug%vecview_solution) then
491) string = 'Surf_xx_' // trim(adjustl(string2)) // '.bin'
492) call PetscViewerBinaryOpen(surf_realization%option%mycomm,string, &
493) FILE_MODE_WRITE,viewer,ierr);CHKERRQ(ierr)
494) call VecView(xx,viewer,ierr);CHKERRQ(ierr)
495) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
496)
497) string = 'Surf_ff_' // trim(adjustl(string2)) // '.bin'
498) call PetscViewerBinaryOpen(surf_realization%option%mycomm,string, &
499) FILE_MODE_WRITE,viewer,ierr);CHKERRQ(ierr)
500) call VecView(ff,viewer,ierr);CHKERRQ(ierr)
501) call PetscViewerDestroy(viewer,ierr);CHKERRQ(ierr)
502)
503) endif
504)
505) end subroutine SurfaceFlowRHSFunction
506)
507) ! ************************************************************************** !
508)
509) subroutine SurfaceFlowComputeMaxDt(surf_realization,max_allowable_dt)
510) !
511) ! This routine maximum allowable 'dt' for surface flow model.
512) !
513) ! Author: Gautam Bisht, LBNL
514) ! Date: 03/07/13
515) !
516)
517)
518) use Connection_module
519) use Realization_Surface_class
520) use Patch_module
521) use Grid_module
522) use Option_module
523) use Coupler_module
524) use Surface_Field_module
525) use Debug_module
526) use Surface_Global_Aux_module
527)
528) implicit none
529)
530) class(realization_surface_type) :: surf_realization
531) PetscErrorCode :: ierr
532)
533) type(grid_type), pointer :: grid
534) type(patch_type), pointer :: patch
535) type(option_type), pointer :: option
536) type(surface_field_type), pointer :: surf_field
537) type(coupler_type), pointer :: boundary_condition
538) type(connection_set_list_type), pointer :: connection_set_list
539) type(connection_set_type), pointer :: cur_connection_set
540) type(surface_global_auxvar_type), pointer :: surf_global_auxvars(:)
541) type(surface_global_auxvar_type), pointer :: surf_global_auxvars_bc(:)
542)
543) PetscInt :: local_id, ghosted_id
544) PetscInt :: local_id_up, local_id_dn
545) PetscInt :: ghosted_id_up, ghosted_id_dn
546) PetscInt :: iconn
547) PetscInt :: sum_connection
548)
549) PetscReal :: dx, dy, dz
550) PetscReal :: dist
551) PetscReal :: vel
552) PetscReal :: slope, slope_dn
553) PetscReal :: rho ! density [kg/m^3]
554) PetscReal :: hw_up, hw_dn ! water height [m]
555) PetscReal :: Res(surf_realization%option%nflowdof), v_darcy
556) PetscReal :: max_allowable_dt
557) PetscReal :: dt
558)
559) PetscReal, pointer :: mannings_loc_p(:),area_p(:)
560) PetscReal, pointer :: xc(:),yc(:),zc(:)
561)
562) patch => surf_realization%patch
563) grid => patch%grid
564) option => surf_realization%option
565) surf_field => surf_realization%surf_field
566) surf_global_auxvars => patch%surf_aux%SurfaceGlobal%auxvars
567) surf_global_auxvars_bc => patch%surf_aux%SurfaceGlobal%auxvars_bc
568)
569) call VecGetArrayF90(surf_field%mannings_loc,mannings_loc_p, &
570) ierr);CHKERRQ(ierr)
571) call VecGetArrayF90(surf_field%area,area_p,ierr);CHKERRQ(ierr)
572)
573) Res = 0.d0
574) max_allowable_dt = 1.d10
575)
576) xc => surf_realization%discretization%grid%x
577) yc => surf_realization%discretization%grid%y
578) zc => surf_realization%discretization%grid%z
579)
580) ! Interior Flux Terms -----------------------------------
581) connection_set_list => grid%internal_connection_set_list
582) cur_connection_set => connection_set_list%first
583) sum_connection = 0
584) do
585) if (.not.associated(cur_connection_set)) exit
586) do iconn = 1, cur_connection_set%num_connections
587) sum_connection = sum_connection + 1
588)
589) ghosted_id_up = cur_connection_set%id_up(iconn)
590) ghosted_id_dn = cur_connection_set%id_dn(iconn)
591)
592) local_id_up = grid%nG2L(ghosted_id_up)
593) local_id_dn = grid%nG2L(ghosted_id_dn)
594)
595) dx = xc(ghosted_id_dn) - xc(ghosted_id_up)
596) dy = yc(ghosted_id_dn) - yc(ghosted_id_up)
597) dz = zc(ghosted_id_dn) - zc(ghosted_id_up)
598) dist = sqrt(dx*dx + dy*dy + dz*dz)
599) slope = dz/dist
600)
601) select case(option%surface_flow_formulation)
602) case (KINEMATIC_WAVE)
603) case (DIFFUSION_WAVE)
604) call SurfaceFlowFlux(surf_global_auxvars(ghosted_id_up), &
605) zc(ghosted_id_up), &
606) mannings_loc_p(ghosted_id_up), &
607) surf_global_auxvars(ghosted_id_dn), &
608) zc(ghosted_id_dn), &
609) mannings_loc_p(ghosted_id_dn), &
610) dist, cur_connection_set%area(iconn), &
611) option,vel,Res)
612) end select
613)
614) patch%internal_velocities(1,sum_connection) = vel
615) patch%internal_flow_fluxes(1,sum_connection) = Res(1)
616) if (abs(vel)>eps) then
617) dt = dist/abs(vel)/4.d0
618) max_allowable_dt = min(max_allowable_dt,dt)
619) endif
620)
621) enddo
622) cur_connection_set => cur_connection_set%next
623) enddo
624)
625) ! Boundary Flux Terms -----------------------------------
626) boundary_condition => patch%boundary_condition_list%first
627) sum_connection = 0
628) do
629) if (.not.associated(boundary_condition)) exit
630)
631) cur_connection_set => boundary_condition%connection_set
632)
633) do iconn = 1, cur_connection_set%num_connections
634) sum_connection = sum_connection + 1
635)
636) local_id = cur_connection_set%id_dn(iconn)
637) ghosted_id_dn = grid%nL2G(local_id)
638)
639) dx = xc(ghosted_id_dn) - cur_connection_set%intercp(1,iconn)
640) dy = yc(ghosted_id_dn) - cur_connection_set%intercp(2,iconn)
641) dz = zc(ghosted_id_dn) - cur_connection_set%intercp(3,iconn)
642) dist = sqrt(dx*dx + dy*dy + dz*dz)
643) slope_dn = dz/dist
644)
645) call SurfaceFlowBCFlux(boundary_condition%flow_condition%itype, &
646) surf_global_auxvars_bc(sum_connection), &
647) slope_dn, &
648) mannings_loc_p(ghosted_id_dn), &
649) cur_connection_set%area(iconn), &
650) option,vel,Res)
651) patch%boundary_velocities(1,sum_connection) = vel
652) patch%boundary_flow_fluxes(1,sum_connection) = Res(1)
653)
654) if (abs(vel)>eps) then
655) dt = dist/abs(vel)/4.d0
656) max_allowable_dt = min(max_allowable_dt,dt)
657) endif
658) enddo
659) boundary_condition => boundary_condition%next
660) enddo
661)
662) call VecRestoreArrayF90(surf_field%mannings_loc,mannings_loc_p, &
663) ierr);CHKERRQ(ierr)
664) call VecRestoreArrayF90(surf_field%area,area_p,ierr);CHKERRQ(ierr)
665)
666) end subroutine SurfaceFlowComputeMaxDt
667)
668) ! ************************************************************************** !
669)
670) subroutine SurfaceFlowFlux(surf_global_auxvar_up, &
671) zc_up, &
672) mannings_up, &
673) surf_global_auxvar_dn, &
674) zc_dn, &
675) mannings_dn, &
676) dist, &
677) length, &
678) option, &
679) vel, &
680) Res)
681) !
682) ! This routine computes the internal flux term for under
683) ! diffusion-wave assumption.
684) !
685) ! Author: Gautam Bisht, LBL
686) ! Date: 08/03/12
687) !
688)
689) use Surface_Global_Aux_module
690) use Option_module
691)
692) implicit none
693)
694) type(option_type) :: option
695) type(surface_global_auxvar_type) :: surf_global_auxvar_up
696) type(surface_global_auxvar_type) :: surf_global_auxvar_dn
697) PetscReal :: zc_up, zc_dn
698) PetscReal :: mannings_up, mannings_dn
699) PetscReal :: head_up, head_dn
700) PetscReal :: dist, length
701) PetscReal :: vel ! units: m/s
702) PetscReal :: Res(1:option%nflowdof) ! units: m^3/s
703)
704) PetscReal :: hw_half
705) PetscReal :: mannings_half
706) PetscReal :: dhead
707)
708) ! upwind the total head and Manning's coefficient
709) head_up = surf_global_auxvar_up%head(1) + zc_up
710) head_dn = surf_global_auxvar_dn%head(1) + zc_dn
711) if (head_up > head_dn) then
712) hw_half = surf_global_auxvar_up%head(1)
713) mannings_half = mannings_up
714) else
715) hw_half = surf_global_auxvar_dn%head(1)
716) mannings_half = mannings_dn
717) endif
718)
719) ! compute Manning's velocity
720) dhead = head_up - head_dn
721) vel = sign(hw_half**(2.d0/3.d0)/mannings_half*abs(dhead/dist)**0.5d0,dhead) ! [m/s]
722)
723) ! compute the volumetric flow rate
724) Res(TH_PRESSURE_DOF) = vel*hw_half*length ! [m^3/s]
725)
726) end subroutine SurfaceFlowFlux
727)
728) ! ************************************************************************** !
729)
730) subroutine SurfaceFlowBCFlux(ibndtype, &
731) surf_global_auxvar, &
732) slope, &
733) mannings, &
734) length, &
735) option, &
736) vel, &
737) Res)
738) !
739) ! This routine computes the boundary term surface water equation
740) !
741) ! Author: Gautam Bisht, LBL
742) ! Date: 08/03/12
743) !
744)
745) use Option_module
746) use Surface_Global_Aux_module
747)
748) implicit none
749)
750) type(option_type) :: option
751) type(surface_global_auxvar_type) :: surf_global_auxvar
752) PetscReal :: slope
753) PetscReal :: mannings
754) PetscReal :: length
755) PetscReal :: flux
756) PetscInt :: ibndtype(:)
757) PetscReal :: vel
758) PetscReal :: Res(1:option%nflowdof)
759)
760) PetscInt :: pressure_bc_type
761) PetscReal :: head
762)
763) flux = 0.d0
764) vel = 0.d0
765)
766) ! Flow
767) pressure_bc_type = ibndtype(TH_PRESSURE_DOF)
768) head = surf_global_auxvar%head(1)
769)
770) select case(pressure_bc_type)
771) case (ZERO_GRADIENT_BC)
772) if (slope<0.d0) then
773) vel = 0.d0
774) else
775) vel = -sqrt(dabs(slope))/mannings*((head)**(2.d0/3.d0))
776) endif
777) case default
778) option%io_buffer = 'Uknown pressure_bc_type for surface flow '
779) end select
780)
781) flux = head*vel
782) Res(RICHARDS_PRESSURE_DOF) = flux*length
783)
784) end subroutine SurfaceFlowBCFlux
785)
786) ! ************************************************************************** !
787)
788) subroutine SurfaceFlowUpdateAuxVars(surf_realization)
789) !
790) ! This routine updates auxiliary variables
791) !
792) ! Author: Gautam Bisht, LBNL
793) ! Date: 03/07/13
794) !
795)
796) use Realization_Surface_class
797) use Patch_module
798) use Option_module
799) use Surface_Field_module
800) use Grid_module
801) use Coupler_module
802) use Connection_module
803) use Surface_Material_module
804) use Surface_Global_Aux_module
805)
806) implicit none
807)
808) class(realization_surface_type) :: surf_realization
809)
810) type(option_type), pointer :: option
811) type(patch_type), pointer :: patch
812) type(grid_type), pointer :: grid
813) type(surface_field_type), pointer :: surf_field
814) type(coupler_type), pointer :: boundary_condition
815) type(coupler_type), pointer :: source_sink
816) type(connection_set_type), pointer :: cur_connection_set
817) type(surface_global_auxvar_type), pointer :: surf_global_auxvars(:)
818) type(surface_global_auxvar_type), pointer :: surf_global_auxvars_bc(:)
819) type(surface_global_auxvar_type), pointer :: surf_global_auxvars_ss(:)
820)
821) PetscInt :: ghosted_id, local_id, sum_connection, iconn
822) PetscReal, pointer :: xx_loc_p(:), icap_loc_p(:), iphase_loc_p(:)
823) PetscReal :: xxbc(1)
824) PetscReal :: xxss(1)
825) PetscReal :: tsrc1
826) PetscErrorCode :: ierr
827)
828) option => surf_realization%option
829) patch => surf_realization%patch
830) grid => patch%grid
831) surf_field => surf_realization%surf_field
832)
833) surf_global_auxvars => patch%surf_aux%SurfaceGlobal%auxvars
834) surf_global_auxvars_bc => patch%surf_aux%SurfaceGlobal%auxvars_bc
835) surf_global_auxvars_ss => patch%surf_aux%SurfaceGlobal%auxvars_ss
836)
837) call VecGetArrayF90(surf_field%flow_xx_loc,xx_loc_p, ierr);CHKERRQ(ierr)
838)
839) ! Internal aux vars
840) do ghosted_id = 1, grid%ngmax
841) if (grid%nG2L(ghosted_id) < 0) cycle ! bypass ghosted corner cells
842)
843) !geh - Ignore inactive cells with inactive materials
844) if (associated(patch%imat)) then
845) if (patch%imat(ghosted_id) <= 0) cycle
846) endif
847) surf_global_auxvars(ghosted_id)%head(1) = xx_loc_p(ghosted_id)
848) enddo
849) call VecRestoreArrayF90(surf_field%flow_xx_loc,xx_loc_p, ierr);CHKERRQ(ierr)
850)
851) call VecGetArrayF90(surf_field%flow_xx_loc,xx_loc_p, ierr);CHKERRQ(ierr)
852) ! Boundary aux vars
853) boundary_condition => patch%boundary_condition_list%first
854) sum_connection = 0
855) do
856) if (.not.associated(boundary_condition)) exit
857) cur_connection_set => boundary_condition%connection_set
858) do iconn = 1, cur_connection_set%num_connections
859) sum_connection = sum_connection + 1
860) local_id = cur_connection_set%id_dn(iconn)
861) ghosted_id = grid%nL2G(local_id)
862) if (associated(patch%imat)) then
863) if (patch%imat(ghosted_id) <= 0) cycle
864) endif
865)
866) select case(boundary_condition%flow_condition%itype(RICHARDS_PRESSURE_DOF))
867) case(DIRICHLET_BC,HYDROSTATIC_BC,SEEPAGE_BC,HET_DIRICHLET)
868) xxbc(1) = boundary_condition%flow_aux_real_var(RICHARDS_PRESSURE_DOF,iconn)
869) case(NEUMANN_BC,ZERO_GRADIENT_BC)
870) xxbc(1) = xx_loc_p(ghosted_id)
871) end select
872)
873) surf_global_auxvars_bc(sum_connection)%head(1) = xxbc(1)
874) enddo
875) boundary_condition => boundary_condition%next
876) enddo
877)
878) ! Source/Sink aux vars
879) ! source/sinks
880) source_sink => patch%source_sink_list%first
881) sum_connection = 0
882) do
883) if (.not.associated(source_sink)) exit
884) cur_connection_set => source_sink%connection_set
885) do iconn = 1, cur_connection_set%num_connections
886) sum_connection = sum_connection + 1
887) local_id = cur_connection_set%id_dn(iconn)
888) ghosted_id = grid%nL2G(local_id)
889) if (patch%imat(ghosted_id) <= 0) cycle
890)
891) xxss = xx_loc_p(ghosted_id)
892) surf_global_auxvars_ss(sum_connection)%head(1) = xxss(1)
893) enddo
894) source_sink => source_sink%next
895) enddo
896)
897) call VecRestoreArrayF90(surf_field%flow_xx_loc,xx_loc_p, ierr);CHKERRQ(ierr)
898)
899) end subroutine SurfaceFlowUpdateAuxVars
900)
901) ! ************************************************************************** !
902)
903) function SurfaceFlowGetTecplotHeader(surf_realization,icolumn)
904) !
905) ! This routine surface flow tecplot file header
906) !
907) ! Author: Gautam Bisht, ORNL
908) ! Date: 05/29/12
909) !
910)
911) use Realization_Surface_class
912) use Option_module
913)
914) implicit none
915)
916) character(len=MAXSTRINGLENGTH) :: SurfaceFlowGetTecplotHeader
917) class(realization_surface_type) :: surf_realization
918) PetscInt :: icolumn
919)
920) character(len=MAXSTRINGLENGTH) :: string, string2
921)
922) string = ''
923)
924) if (icolumn > -1) then
925) icolumn = icolumn + 1
926) write(string2,'('',"'',i2,''-P [Pa]"'')') icolumn
927) else
928) write(string2,'('',"P [m]"'')')
929) endif
930) string = trim(string) // trim(string2)
931) #ifdef GLENN_NEW_IO
932) !call OutputOptionAddPlotVariable(realization%output_option,PRESSURE, &
933) ! ZERO_INTEGER,ZERO_INTEGER)
934) #endif
935)
936) SurfaceFlowGetTecplotHeader = string
937)
938) end function SurfaceFlowGetTecplotHeader
939)
940) ! ************************************************************************** !
941)
942) subroutine SurfaceFlowUpdateSurfState(surf_realization)
943) !
944) ! This routine gets updated values of standing water at the end of
945) ! subsurface-flow model timestep.
946) !
947) ! Author: Gautam Bisht, LBL
948) ! Date: 07/30/13
949) !
950)
951) use Connection_module
952) use Coupler_module
953) use Discretization_module
954) use DM_Kludge_module
955) use Grid_module
956) use Option_module
957) use Patch_module
958) use Realization_Subsurface_class
959) use Realization_Base_class
960) use String_module
961) use Surface_Field_module
962) use Realization_Surface_class
963) use EOS_Water_module
964)
965) implicit none
966)
967) #include "petsc/finclude/petscvec.h"
968) #include "petsc/finclude/petscvec.h90"
969) #include "petsc/finclude/petscmat.h"
970) #include "petsc/finclude/petscmat.h90"
971)
972) class(realization_surface_type) :: surf_realization
973)
974) type(coupler_list_type), pointer :: coupler_list
975) type(coupler_type), pointer :: coupler
976) type(connection_set_type), pointer :: cur_connection_set
977) type(dm_ptr_type), pointer :: dm_ptr
978) type(grid_type),pointer :: grid,surf_grid
979) type(option_type), pointer :: option
980) type(patch_type),pointer :: patch,surf_patch
981) type(surface_field_type),pointer :: surf_field
982)
983) PetscInt :: iconn
984) PetscInt :: local_id
985) PetscInt :: sum_connection
986)
987) PetscReal :: den
988) PetscReal :: dum1
989) PetscReal, pointer :: avg_vdarcy_p(:) ! avg darcy velocity [m/s]
990) PetscReal, pointer :: hw_p(:) ! head [m]
991) PetscReal, pointer :: surfpress_p(:)
992) PetscErrorCode :: ierr
993)
994) PetscBool :: coupler_found = PETSC_FALSE
995)
996) option => surf_realization%option
997) surf_field => surf_realization%surf_field
998) surf_grid => surf_realization%discretization%grid
999)
1000) call EOSWaterdensity(option%reference_temperature, &
1001) option%reference_pressure,den,dum1,ierr)
1002)
1003) call VecGetArrayF90(surf_field%flow_xx, hw_p, ierr);CHKERRQ(ierr)
1004) call VecGetArrayF90(surf_field%press_subsurf, surfpress_p, &
1005) ierr);CHKERRQ(ierr)
1006)
1007) do local_id = 1,surf_grid%nlmax
1008)
1009) hw_p(local_id) = (surfpress_p(local_id)-option%reference_pressure)/ &
1010) (abs(option%gravity(3)))/den
1011) if (hw_p(local_id)<1.d-15) hw_p(local_id) = 0.d0
1012)
1013) enddo
1014) call VecRestoreArrayF90(surf_field%flow_xx, hw_p, ierr);CHKERRQ(ierr)
1015) call VecRestoreArrayF90(surf_field%press_subsurf, surfpress_p, &
1016) ierr);CHKERRQ(ierr)
1017)
1018) call DiscretizationGlobalToLocal(surf_realization%discretization, &
1019) surf_field%flow_xx, &
1020) surf_field%flow_xx_loc, &
1021) NFLOWDOF)
1022) call SurfaceFlowUpdateAuxVars(surf_realization)
1023)
1024) end subroutine SurfaceFlowUpdateSurfState
1025)
1026) ! ************************************************************************** !
1027)
1028) end module Surface_Flow_module