Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
direct_net_infiltration__gridded_data.F90
Go to the documentation of this file.
1!> @file
2!! Contains the module \ref direct_net_infiltration__gridded_data.
3
4!>
5!! Module \ref direct_net_infiltration__gridded_data
6!! provides support for adding miscellaneous source and sink terms.
7
9
10 use iso_c_binding
12 use data_catalog
14 use datetime
15 use dictionary
16 use exceptions
19 use parameters, only : params
21 use fstring
22 use fstring_list
23
24 implicit none
25
26 private
27
29
35
36 real (c_float), allocatable :: fcesspool(:)
37 real (c_float), allocatable :: fdisposal_well(:)
38 real (c_float), allocatable :: fwater_body_recharge(:)
39 real (c_float), allocatable :: fwater_main(:)
40 real (c_float), allocatable :: fannual_recharge_rate(:)
41
42 ! ****_TABLE variables: will have same number of values as there are landuses
43 real (c_float), allocatable :: fcesspool_table(:)
44 real (c_float), allocatable :: fdisposal_well_table(:)
45 real (c_float), allocatable :: fwater_body_recharge_table(:)
46 real (c_float), allocatable :: fwater_main_table(:)
47 real (c_float), allocatable :: fannual_recharge_rate_table(:)
48
49 type (t_netcdf4_file), pointer :: pncfile
50
52
53contains
54
55 !> Initialize the routine to enable input/output of arbitrary sources/sink terms.
56 !!
57 !! Open gridded data file.
58 !! Open a NetCDF output file to hold variable output.
59 !!
60 !! @param[in] is_cell_active 2D array of active cells within the model domain.
61 !! @param[in] iLanduseIndex 1D vector of indices corresponding to rows of the
62 !! landuse lookup table(s).
63 !! @param[in] dX 1D vector of X coordinates associated with the model domain.
64 !! @param[in] dY 1D vector of Y coordinates.
65 !! @param[in] dX_lon 2D array of longitude values.
66 !! @param[in] dY_lat 2D array of latitude values.
67
68 subroutine direct_net_infiltration_initialize( is_cell_active, landuse_index )
69
70 logical (c_bool), intent(in) :: is_cell_active(:,:)
71 integer (c_int), intent(in) :: landuse_index(:)
72
73 ! [ LOCALS ]
74 integer (c_int) :: status
75 type (fstring_list_t) :: parameter_list
76 integer (c_int) :: indx
77 integer (c_int) :: inx
78 integer (c_int) :: iny
79 integer (c_int), allocatable :: landuse_codes(:)
80 integer (c_int) :: number_of_landuses
81 logical (c_bool) :: are_lengths_equal
82
83
84 !> Determine how many landuse codes are present
85 call parameter_list%append( "LU_Code" )
86 call parameter_list%append( "Landuse_Code" )
87 call parameter_list%append( "Landuse_Lookup_Code" )
88
89 call params%get_parameters( slkeys=parameter_list, ivalues=landuse_codes )
90 number_of_landuses = count( landuse_codes > 0 )
91
92 call parameter_list%clear()
93 call parameter_list%append( "Annual_direct_net_infiltration_rate" )
94 call parameter_list%append( "Annual_direct_recharge_rate" )
95 call parameter_list%append( "Annual_direct_recharge" )
96 call parameter_list%append( "Annual_direct_net_infiltration" )
97
98 call params%get_parameters( slkeys=parameter_list , fvalues=fannual_recharge_rate_table )
99
100 ! attempt to find a source of GRIDDED ANNUAL DIRECT RECHARGE data
101 pannual_recharge_rate => dat%find( "ANNUAL_DIRECT_NET_INFILTRATION_RATE" )
102
103
104 call parameter_list%clear()
105 call parameter_list%append( "Cesspool_direct_net_infiltration" )
106 call parameter_list%append( "Cesspool_recharge" )
107 call parameter_list%append( "Cesspool_discharge" )
108 call parameter_list%append( "Cesspool_leakage" )
109
110 call params%get_parameters( slkeys=parameter_list , fvalues=fcesspool_table )
111
112 ! attempt to find a source of GRIDDED CESSPOOL data
113 pcesspool => dat%find( "CESSPOOL_LEAKAGE" )
114
115
116 ! call parameter_list%clear()
117 ! call parameter_list%append( "Storm_drain_discharge" )
118 ! call parameter_list%append( "Storm_drain_recharge" )
119 ! call parameter_list%append( "Storm_drain_leakage" )
120 !
121 ! call PARAMS%get_parameters( slKeys=parameter_list , fValues=fSTORM_DRAIN_TABLE )
122 !
123 ! pSTORM_DRAIN => DAT%find( "STORM_DRAIN" )
124 !
125
126 call parameter_list%clear()
127 call parameter_list%append( "Water_body_recharge" )
128 call parameter_list%append( "Water_body_discharge" )
129 call parameter_list%append( "Water_body_leakage" )
130
131 call params%get_parameters( slkeys=parameter_list , fvalues=fwater_body_recharge_table )
132
133 pwater_body_recharge => dat%find( "WATER_BODY_RECHARGE" )
134
135
136 call parameter_list%clear()
137 call parameter_list%append( "Water_main_recharge" )
138 call parameter_list%append( "Water_main_discharge" )
139 call parameter_list%append( "Water_main_leakage" )
140
141 call params%get_parameters( slkeys=parameter_list , fvalues=fwater_main_table )
142
143 pwater_main => dat%find( "WATER_MAIN_LEAKAGE" )
144
145
146 call parameter_list%clear()
147 call parameter_list%append( "Disposal_well_recharge" )
148 call parameter_list%append( "Disposal_well_discharge" )
149
150 call params%get_parameters( slkeys=parameter_list , fvalues=fdisposal_well_table )
151
152 pdisposal_well => dat%find( "DISPOSAL_WELL_DISCHARGE" )
153
154 if ( associated( pannual_recharge_rate ) ) then
155
156 allocate( fannual_recharge_rate( count( is_cell_active ) ), stat=status )
157 call assert( status==0, "Problem allocating memory", __file__, __line__ )
158
159 elseif ( fannual_recharge_rate_table(1) > ftinyval ) then
160
161 are_lengths_equal = ( ( ubound(fannual_recharge_rate_table,1) == ubound(landuse_codes,1) ) )
162
163 if ( .not. are_lengths_equal ) &
164 call warn( smessage="The number of landuses does not match the number of annual direct" &
165 //" recharge rate values.", smodule=__file__, iline=__line__, lfatal=.true._c_bool )
166
167 allocate( fannual_recharge_rate( count( is_cell_active ) ), stat=status )
168 call assert( status==0, "Problem allocating memory", __file__, __line__ )
169
170 do indx=lbound( landuse_index, 1 ), ubound( landuse_index, 1 )
171 fannual_recharge_rate( indx ) = fannual_recharge_rate_table( landuse_index( indx ) )
172 enddo
173
174 endif
175
176
177 if ( associated( pcesspool ) ) then
178
179 allocate( fcesspool( count( is_cell_active ) ), stat=status )
180 call assert( status==0, "Problem allocating memory", __file__, __line__ )
181
182 elseif ( fcesspool_table(1) > ftinyval ) then
183
184 are_lengths_equal = ( ( ubound(fcesspool_table,1) == ubound(landuse_codes,1) ) )
185
186 if ( .not. are_lengths_equal ) &
187 call warn( smessage="The number of landuses does not match the number of cesspool discharge/leakage values.", &
188 smodule=__file__, iline=__line__, lfatal=.true._c_bool )
189
190 allocate( fcesspool( count( is_cell_active ) ), stat=status )
191 call assert( status==0, "Problem allocating memory", __file__, __line__ )
192
193 do indx=lbound( landuse_index, 1 ), ubound( landuse_index, 1 )
194 fcesspool( indx ) = fcesspool_table( landuse_index( indx ) )
195 enddo
196
197 endif
198
199 ! if ( associated( pSTORM_DRAIN ) ) then
200 !
201 ! allocate( fSTORM_DRAIN( count( is_cell_active ) ), stat=status )
202 ! call assert( status==0, "Problem allocating memory", __FILE__, __LINE__ )
203 !
204 ! elseif ( fSTORM_DRAIN_TABLE(1) > fTINYVAL ) then
205 !
206 ! are_lengths_equal = ( ( ubound(fSTORM_DRAIN_TABLE,1) == ubound(landuse_codes,1) ) )
207 !
208 ! if ( .not. are_lengths_equal ) &
209 ! call warn( sMessage="The number of landuses does not match the number of storm drain discharge/leakage values.", &
210 ! sModule=__FILE__, iLine=__LINE__, lFatal=.true._c_bool )
211 !
212 ! allocate( fSTORM_DRAIN( count( is_cell_active ) ), stat=status )
213 ! call assert( status==0, "Problem allocating memory", __FILE__, __LINE__ )
214 !
215 ! do indx=lbound( landuse_index, 1 ), ubound( landuse_index, 1 )
216 ! fSTORM_DRAIN( indx ) = fSTORM_DRAIN_TABLE( landuse_index( indx ) )
217 ! enddo
218 !
219 ! endif
220
221 if ( associated( pwater_body_recharge ) ) then
222
223 allocate( fwater_body_recharge( count( is_cell_active ) ), stat=status )
224 call assert( status==0, "Problem allocating memory", __file__, __line__ )
225
226 elseif ( fwater_body_recharge_table(1) > ftinyval ) then
227
228 are_lengths_equal = ( ( ubound(fwater_body_recharge_table,1) == ubound(landuse_codes,1) ) )
229
230 if ( .not. are_lengths_equal ) &
231 call warn( smessage="The number of landuses does not match the number of water body recharge/leakage values.", &
232 smodule=__file__, iline=__line__, lfatal=.true._c_bool )
233
234 allocate( fwater_body_recharge( count( is_cell_active ) ), stat=status )
235 call assert( status==0, "Problem allocating memory", __file__, __line__ )
236
237 do indx=lbound( landuse_index, 1 ), ubound( landuse_index, 1 )
238 fwater_body_recharge( indx ) = fwater_body_recharge_table( landuse_index( indx ) )
239 enddo
240
241 endif
242
243
244 if ( associated( pwater_main ) ) then
245
246 allocate( fwater_main( count( is_cell_active ) ), stat=status )
247 call assert( status==0, "Problem allocating memory", __file__, __line__ )
248
249 elseif ( fwater_main_table(1) > ftinyval ) then
250
251 are_lengths_equal = ( ( ubound(fwater_main_table,1) == ubound(landuse_codes,1) ) )
252
253 if ( .not. are_lengths_equal ) &
254 call warn( smessage="The number of landuses does not match the number of water main leakage values.", &
255 smodule=__file__, iline=__line__, lfatal=.true._c_bool )
256
257 allocate( fwater_main( count( is_cell_active ) ), stat=status )
258 call assert( status==0, "Problem allocating memory", __file__, __line__ )
259
260 do indx=lbound( landuse_index, 1 ), ubound( landuse_index, 1 )
261 fwater_main( indx ) = fwater_main_table( landuse_index( indx ) )
262 enddo
263
264 endif
265
266
267 if ( associated( pdisposal_well ) ) then
268
269 allocate( fdisposal_well( count( is_cell_active ) ), stat=status )
270 call assert( status==0, "Problem allocating memory", __file__, __line__ )
271
272 elseif ( fdisposal_well_table(1) > ftinyval ) then
273
274 are_lengths_equal = ( ( ubound(fdisposal_well_table,1) == ubound(landuse_codes,1) ) )
275
276 if ( .not. are_lengths_equal ) &
277 call warn( smessage="The number of landuses does not match the number of discharge well values.", &
278 smodule=__file__, iline=__line__, lfatal=.true._c_bool )
279
280 allocate( fdisposal_well( count( is_cell_active ) ), stat=status )
281 call assert( status==0, "Problem allocating memory", __file__, __line__ )
282
283 do indx=lbound( landuse_index, 1 ), ubound( landuse_index, 1 )
284 fdisposal_well( indx ) = fdisposal_well_table( landuse_index( indx ) )
285 enddo
286
287 endif
288
289 ! initialize last retrieval date to something implausibly low to trigger initial read
290 ! in the calculate procedure
291 call date_of_last_retrieval%parseDate("01/01/1000", __file__, __line__)
292
294
295!--------------------------------------------------------------------------------------------------
296
297 subroutine direct_net_infiltration_calculate( direct_net_infiltration, indx, is_cell_active, nodata_fill_value )
298
299 real (c_float), intent(inout) :: direct_net_infiltration
300 integer (c_int), intent(in) :: indx
301 logical (c_bool), intent(in) :: is_cell_active(:,:)
302 real (c_float), intent(in) :: nodata_fill_value(:,:)
303
304 ! [ LOCALS ]
305 real (c_float) :: ffactor
306
307 if ( .not. date_of_last_retrieval == sim_dt%curr ) then
308
309 associate( dt => sim_dt%curr )
310
311 if ( associated( pcesspool ) ) then
312 call pcesspool%getvalues( dt )
313 if ( pcesspool%lGridHasChanged ) fcesspool = pack( pcesspool%pGrdBase%rData, is_cell_active )
314 endif
315
316 if ( associated( pdisposal_well ) ) then
317 call pdisposal_well%getvalues( dt )
318 if ( pdisposal_well%lGridHasChanged ) fdisposal_well = &
319 pack( pdisposal_well%pGrdBase%rData, is_cell_active )
320 endif
321
322 ! if ( associated( pSTORM_DRAIN ) ) then
323 ! call pSTORM_DRAIN%getvalues( iMonth, iDay, iYear, iJulianDay )
324 ! if ( pSTORM_DRAIN%lGridHasChanged ) fSTORM_DRAIN = pack( pSTORM_DRAIN%pGrdBase%rData, is_cell_active )
325 ! endif
326
327 if ( associated( pwater_body_recharge ) ) then
328 call pwater_body_recharge%getvalues( dt )
329 if ( pwater_body_recharge%lGridHasChanged ) fwater_body_recharge = &
330 pack( pwater_body_recharge%pGrdBase%rData, is_cell_active )
331 endif
332
333 if ( associated( pwater_main ) ) then
334 call pwater_main%getvalues( dt )
335 if ( pwater_main%lGridHasChanged ) fwater_main = pack( pwater_main%pGrdBase%rData, is_cell_active )
336 endif
337
338 if ( associated( pannual_recharge_rate ) ) then
339 call pannual_recharge_rate%getvalues( dt )
340 if ( pannual_recharge_rate%lGridHasChanged ) fannual_recharge_rate = &
341 pack( pannual_recharge_rate%pGrdBase%rData, is_cell_active )
342 endif
343
345
346 end associate
347
348 endif
349
350 direct_net_infiltration = 0.0_c_float
351
352 if ( allocated( fcesspool ) ) &
353 direct_net_infiltration = direct_net_infiltration + fcesspool( indx )
354
355 if ( allocated( fdisposal_well ) ) &
356 direct_net_infiltration = direct_net_infiltration + fdisposal_well( indx )
357
358 if ( allocated( fwater_main ) ) &
359 direct_net_infiltration = direct_net_infiltration + fwater_main( indx )
360
361 if ( allocated( fwater_body_recharge ) ) &
362 direct_net_infiltration = direct_net_infiltration + fwater_body_recharge( indx )
363
364! if ( allocated( fSTORM_DRAIN ) ) direct_net_infiltration = direct_net_infiltration + fSTORM_DRAIN( indx )
365
366 if ( allocated( fannual_recharge_rate) ) &
367 direct_net_infiltration = direct_net_infiltration + fannual_recharge_rate( indx ) / 365.25_c_float
368
370
This module contains physical constants and convenience functions aimed at performing unit conversion...
logical(c_bool), parameter, public true
Defines the DATA_CATALOG_T data type, which contains type-bound procedures to add,...
type(data_catalog_t), public dat
DAT is a global to hold data catalog entries.
This module contains the DATETIME_T class and associated time and date-related routines,...
Definition datetime.F90:9
Module direct_net_infiltration__gridded_data provides support for adding miscellaneous source and sin...
real(c_float), dimension(:), allocatable fwater_body_recharge
real(c_float), dimension(:), allocatable fdisposal_well_table
real(c_float), dimension(:), allocatable fdisposal_well
real(c_float), dimension(:), allocatable fwater_main_table
real(c_float), dimension(:), allocatable fwater_body_recharge_table
subroutine, public direct_net_infiltration_calculate(direct_net_infiltration, indx, is_cell_active, nodata_fill_value)
real(c_float), dimension(:), allocatable fannual_recharge_rate_table
real(c_float), dimension(:), allocatable fwater_main
real(c_float), dimension(:), allocatable fcesspool_table
real(c_float), dimension(:), allocatable fannual_recharge_rate
subroutine, public direct_net_infiltration_initialize(is_cell_active, landuse_index)
Initialize the routine to enable input/output of arbitrary sources/sink terms.
subroutine, public warn(smessage, smodule, iline, shints, lfatal, iloglevel, lecho)
Provide support for use of netCDF files as input for time-varying, gridded meteorlogic data,...
type(parameters_t), public params
type(date_range_t), public sim_dt