MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
virtualbasemodule Module Reference

Data Types

type  virtualdatatype
 This is a generic data structure to virtualize pieces of memory in 2 distinct ways: More...
 
type  virtualinttype
 
type  virtualint1dtype
 
type  virtualdbltype
 
type  virtualdbl1dtype
 
type  virtualdbl2dtype
 
interface  vm_allocate_if
 
interface  vm_deallocate_if
 

Functions/Subroutines

class(virtualdatatype) function, pointer vm_to_base (this)
 
logical(lgp) function vm_check_stage (this, stage)
 Check if this data item requires syncing. More...
 
subroutine vm_link (this)
 
integer(i4b) function, dimension(:), pointer, contiguous get_element_map (this)
 Return array with offsets for elements. More...
 
subroutine vm_allocate_int (this, var_name, mem_path, shape)
 
subroutine vm_deallocate_int (this)
 
subroutine vm_allocate_int1d (this, var_name, mem_path, shape)
 
subroutine vm_deallocate_int1d (this)
 
subroutine vm_allocate_dbl (this, var_name, mem_path, shape)
 
subroutine vm_deallocate_dbl (this)
 
subroutine vm_allocate_dbl1d (this, var_name, mem_path, shape)
 
subroutine vm_deallocate_dbl1d (this)
 
subroutine vm_allocate_dbl2d (this, var_name, mem_path, shape)
 
subroutine vm_deallocate_dbl2d (this)
 
integer(i4b) function get_int (this)
 
integer(i4b) function get_int1d (this, i_rmt)
 
integer(i4b) function, dimension(:), pointer, contiguous get_array_int1d (this)
 
real(dp) function get_dbl (this)
 
real(dp) function get_dbl1d (this, i_rmt)
 
real(dp) function, dimension(:), pointer, contiguous get_array_dbl1d (this)
 
real(dp) function get_dbl2d (this, j_cmp, i_rmt)
 
real(dp) function, dimension(:, :), pointer, contiguous get_array_dbl2d (this)
 
class(virtualdatatype) function, pointer, public get_virtual_data_from_list (list, idx)
 

Variables

integer(i4b), parameter, public map_all_type = 0
 
integer(i4b), parameter, public map_node_type = 1
 
integer(i4b), parameter, public map_conn_type = 2
 
integer(i4b), parameter, public nr_vdc_element_maps = 2
 

Function/Subroutine Documentation

◆ get_array_dbl1d()

real(dp) function, dimension(:), pointer, contiguous virtualbasemodule::get_array_dbl1d ( class(virtualdbl1dtype this)
private

Definition at line 307 of file VirtualBase.f90.

308  class(VirtualDbl1dType) :: this
309  real(DP), dimension(:), pointer, contiguous :: array
310 
311  array => this%virtual_mt%adbl1d
312 

◆ get_array_dbl2d()

real(dp) function, dimension(:, :), pointer, contiguous virtualbasemodule::get_array_dbl2d ( class(virtualdbl2dtype this)
private

Definition at line 332 of file VirtualBase.f90.

333  class(VirtualDbl2dType) :: this
334  real(DP), dimension(:, :), pointer, contiguous :: array
335 
336  array => this%virtual_mt%adbl2d
337 

◆ get_array_int1d()

integer(i4b) function, dimension(:), pointer, contiguous virtualbasemodule::get_array_int1d ( class(virtualint1dtype this)
private

Definition at line 275 of file VirtualBase.f90.

276  class(VirtualInt1dType) :: this
277  integer(I4B), dimension(:), pointer, contiguous :: array
278 
279  array => this%virtual_mt%aint1d
280 

◆ get_dbl()

real(dp) function virtualbasemodule::get_dbl ( class(virtualdbltype this)
private

Definition at line 283 of file VirtualBase.f90.

284  class(VirtualDblType) :: this
285  real(DP) :: val
286 
287  val = this%virtual_mt%dblsclr
288 

◆ get_dbl1d()

real(dp) function virtualbasemodule::get_dbl1d ( class(virtualdbl1dtype this,
integer(i4b)  i_rmt 
)
private

Definition at line 291 of file VirtualBase.f90.

292  class(VirtualDbl1dType) :: this
293  integer(I4B) :: i_rmt
294  real(DP) :: val
295  ! local
296  integer(I4B) :: i_vrt
297 
298  if (this%is_reduced) then
299  i_vrt = this%remote_to_virtual(i_rmt)
300  else
301  i_vrt = i_rmt
302  end if
303  val = this%virtual_mt%adbl1d(i_vrt)
304 

◆ get_dbl2d()

real(dp) function virtualbasemodule::get_dbl2d ( class(virtualdbl2dtype this,
integer(i4b)  j_cmp,
integer(i4b)  i_rmt 
)
private

Definition at line 315 of file VirtualBase.f90.

316  class(VirtualDbl2dType) :: this
317  integer(I4B) :: j_cmp
318  integer(I4B) :: i_rmt
319  real(DP) :: val
320  ! local
321  integer(I4B) :: i_vrt
322 
323  if (this%is_reduced) then
324  i_vrt = this%remote_to_virtual(i_rmt)
325  else
326  i_vrt = i_rmt
327  end if
328  val = this%virtual_mt%adbl2d(j_cmp, i_vrt)
329 

◆ get_element_map()

integer(i4b) function, dimension(:), pointer, contiguous virtualbasemodule::get_element_map ( class(virtualdatatype), target  this)
private

Definition at line 155 of file VirtualBase.f90.

156  class(VirtualDataType), target :: this
157  integer(I4B), dimension(:), pointer, contiguous :: el_map
158 
159  el_map => null()
160  if (this%map_type > 0) then
161  el_map => this%remote_elem_shift
162  end if
163 

◆ get_int()

integer(i4b) function virtualbasemodule::get_int ( class(virtualinttype this)
private

Definition at line 251 of file VirtualBase.f90.

252  class(VirtualIntType) :: this
253  integer(I4B) :: val
254 
255  val = this%virtual_mt%intsclr
256 

◆ get_int1d()

integer(i4b) function virtualbasemodule::get_int1d ( class(virtualint1dtype this,
integer(i4b)  i_rmt 
)
private

Definition at line 259 of file VirtualBase.f90.

260  class(VirtualInt1dType) :: this
261  integer(I4B) :: i_rmt
262  integer(I4B) :: val
263  ! local
264  integer(I4B) :: i_vrt
265 
266  if (this%is_reduced) then
267  i_vrt = this%remote_to_virtual(i_rmt)
268  else
269  i_vrt = i_rmt
270  end if
271  val = this%virtual_mt%aint1d(i_vrt)
272 

◆ get_virtual_data_from_list()

class(virtualdatatype) function, pointer, public virtualbasemodule::get_virtual_data_from_list ( type(listtype list,
integer(i4b)  idx 
)

Definition at line 340 of file VirtualBase.f90.

341  type(ListType) :: list
342  integer(I4B) :: idx
343  class(VirtualDataType), pointer :: vd
344  ! local
345  class(*), pointer :: obj_ptr
346 
347  vd => null()
348  obj_ptr => list%GetItem(idx)
349  select type (obj_ptr)
350  class is (virtualdatatype)
351  vd => obj_ptr
352  end select
353 
Here is the caller graph for this function:

◆ vm_allocate_dbl()

subroutine virtualbasemodule::vm_allocate_dbl ( class(virtualdbltype this,
character(len=*)  var_name,
character(len=*)  mem_path,
integer(i4b), dimension(:)  shape 
)
private

Definition at line 200 of file VirtualBase.f90.

201  class(VirtualDblType) :: this
202  character(len=*) :: var_name
203  character(len=*) :: mem_path
204  integer(I4B), dimension(:) :: shape
205 
206  call mem_allocate(this%dblsclr, var_name, mem_path)
207 

◆ vm_allocate_dbl1d()

subroutine virtualbasemodule::vm_allocate_dbl1d ( class(virtualdbl1dtype this,
character(len=*)  var_name,
character(len=*)  mem_path,
integer(i4b), dimension(:)  shape 
)
private

Definition at line 217 of file VirtualBase.f90.

218  class(VirtualDbl1dType) :: this
219  character(len=*) :: var_name
220  character(len=*) :: mem_path
221  integer(I4B), dimension(:) :: shape
222 
223  call mem_allocate(this%dbl1d, shape(1), var_name, mem_path)
224 

◆ vm_allocate_dbl2d()

subroutine virtualbasemodule::vm_allocate_dbl2d ( class(virtualdbl2dtype this,
character(len=*)  var_name,
character(len=*)  mem_path,
integer(i4b), dimension(:)  shape 
)
private

Definition at line 234 of file VirtualBase.f90.

235  class(VirtualDbl2dType) :: this
236  character(len=*) :: var_name
237  character(len=*) :: mem_path
238  integer(I4B), dimension(:) :: shape
239 
240  call mem_allocate(this%dbl2d, shape(1), shape(2), var_name, mem_path)
241 

◆ vm_allocate_int()

subroutine virtualbasemodule::vm_allocate_int ( class(virtualinttype this,
character(len=*)  var_name,
character(len=*)  mem_path,
integer(i4b), dimension(:)  shape 
)
private

Definition at line 166 of file VirtualBase.f90.

167  class(VirtualIntType) :: this
168  character(len=*) :: var_name
169  character(len=*) :: mem_path
170  integer(I4B), dimension(:) :: shape
171 
172  call mem_allocate(this%intsclr, var_name, mem_path)
173 

◆ vm_allocate_int1d()

subroutine virtualbasemodule::vm_allocate_int1d ( class(virtualint1dtype this,
character(len=*)  var_name,
character(len=*)  mem_path,
integer(i4b), dimension(:)  shape 
)
private

Definition at line 183 of file VirtualBase.f90.

184  class(VirtualInt1dType) :: this
185  character(len=*) :: var_name
186  character(len=*) :: mem_path
187  integer(I4B), dimension(:) :: shape
188 
189  call mem_allocate(this%int1d, shape(1), var_name, mem_path)
190 

◆ vm_check_stage()

logical(lgp) function virtualbasemodule::vm_check_stage ( class(virtualdatatype), target  this,
integer(i4b)  stage 
)
private

Definition at line 129 of file VirtualBase.f90.

130  use arrayhandlersmodule, only: ifind
131  class(VirtualDataType), target :: this
132  integer(I4B) :: stage, stg_idx
133  logical(LGP) :: has_stage
134 
135  has_stage = .false.
136  if (allocated(this%sync_stages)) then
137  stg_idx = ifind(this%sync_stages, stage)
138  has_stage = (stg_idx > 0)
139  end if
140 

◆ vm_deallocate_dbl()

subroutine virtualbasemodule::vm_deallocate_dbl ( class(virtualdbltype this)
private

Definition at line 210 of file VirtualBase.f90.

211  class(VirtualDblType) :: this
212 
213  if (this%is_remote) call mem_deallocate(this%dblsclr)
214 

◆ vm_deallocate_dbl1d()

subroutine virtualbasemodule::vm_deallocate_dbl1d ( class(virtualdbl1dtype this)
private

Definition at line 227 of file VirtualBase.f90.

228  class(VirtualDbl1dType) :: this
229 
230  if (this%is_remote) call mem_deallocate(this%dbl1d)
231 

◆ vm_deallocate_dbl2d()

subroutine virtualbasemodule::vm_deallocate_dbl2d ( class(virtualdbl2dtype this)
private

Definition at line 244 of file VirtualBase.f90.

245  class(VirtualDbl2dType) :: this
246 
247  if (this%is_remote) call mem_deallocate(this%dbl2d)
248 

◆ vm_deallocate_int()

subroutine virtualbasemodule::vm_deallocate_int ( class(virtualinttype this)
private

Definition at line 176 of file VirtualBase.f90.

177  class(VirtualIntType) :: this
178 
179  if (this%is_remote) call mem_deallocate(this%intsclr)
180 

◆ vm_deallocate_int1d()

subroutine virtualbasemodule::vm_deallocate_int1d ( class(virtualint1dtype this)
private

Definition at line 193 of file VirtualBase.f90.

194  class(VirtualInt1dType) :: this
195 
196  if (this%is_remote) call mem_deallocate(this%int1d)
197 

◆ vm_link()

subroutine virtualbasemodule::vm_link ( class(virtualdatatype), target  this)

Definition at line 143 of file VirtualBase.f90.

144  class(VirtualDataType), target :: this
145  ! local
146  logical(LGP) :: found
147 
148  call get_from_memorystore(this%var_name, this%mem_path, &
149  this%virtual_mt, found)
150 
Here is the call graph for this function:

◆ vm_to_base()

class(virtualdatatype) function, pointer virtualbasemodule::vm_to_base ( class(virtualdatatype), target  this)
private

Definition at line 119 of file VirtualBase.f90.

120  class(VirtualDataType), target :: this
121  class(VirtualDataType), pointer :: base_ptr
122 
123  base_ptr => this
124 

Variable Documentation

◆ map_all_type

integer(i4b), parameter, public virtualbasemodule::map_all_type = 0

Definition at line 13 of file VirtualBase.f90.

13  integer(I4B), public, parameter :: MAP_ALL_TYPE = 0

◆ map_conn_type

integer(i4b), parameter, public virtualbasemodule::map_conn_type = 2

Definition at line 15 of file VirtualBase.f90.

15  integer(I4B), public, parameter :: MAP_CONN_TYPE = 2

◆ map_node_type

integer(i4b), parameter, public virtualbasemodule::map_node_type = 1

Definition at line 14 of file VirtualBase.f90.

14  integer(I4B), public, parameter :: MAP_NODE_TYPE = 1

◆ nr_vdc_element_maps

integer(i4b), parameter, public virtualbasemodule::nr_vdc_element_maps = 2

Definition at line 16 of file VirtualBase.f90.

16  integer(I4B), public, parameter :: NR_VDC_ELEMENT_MAPS = 2