MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
VirtualDataContainer.f90
Go to the documentation of this file.
3  use simmodule, only: ustop
4  use listmodule
5  use kindmodule, only: i4b, lgp
10  implicit none
11  private
12 
13  public :: get_vdc_from_list
14  public :: vdc_type_to_str
15 
16  integer(I4B), public, parameter :: vdc_unknown_type = 0
17  integer(I4B), public, parameter :: vdc_gwfmodel_type = 1
18  integer(I4B), public, parameter :: vdc_gwtmodel_type = 2
19  integer(I4B), public, parameter :: vdc_gwemodel_type = 3
20  integer(I4B), public, parameter :: vdc_gwfexg_type = 4
21  integer(I4B), public, parameter :: vdc_gwtexg_type = 5
22  integer(I4B), public, parameter :: vdc_gweexg_type = 6
23  integer(I4B), public, parameter :: vdc_gwfmvr_type = 7
24  integer(I4B), public, parameter :: vdc_gwtmvt_type = 8
25  integer(I4B), public, parameter :: vdc_gwemve_type = 9
26 
27  !> @brief Wrapper for virtual data containers
28  !!
29  !! We can't have an array of pointers in Fortran, so we use
30  !! this trick where we wrap the pointer and have an array
31  !< of VdcPtrType instead.
32  type, public :: vdcptrtype
33  class(virtualdatacontainertype), pointer :: ptr => null()
34  end type vdcptrtype
35 
36  type, public :: vdcelementmaptype
37  integer(I4B) :: nr_virt_elems !< nr. of virtualized elements
38  integer(I4B), dimension(:), pointer, contiguous :: remote_elem_shift => null() !< array with 0-based remote indexes
39  end type vdcelementmaptype
40 
42  integer(I4B) :: max_remote_idx !< max. remote index, also size of the lookup table
43  integer(I4B), dimension(:), pointer, contiguous :: remote_to_virtual => null() !< (sparse) array with local indexes
44  end type vdcelementluttype
45 
46  !> @brief Container (list) of virtual data items.
47  !!
48  !! A virtual model or exchange derives from this base
49  !! and can add the component-specific items to the list
50  !! of virtual data items. As far as synchronization
51  !! of virtual objects is concerned, all that is needed
52  !< is the list of virtual data items in this container.
53  type, public :: virtualdatacontainertype
54  integer(I4B) :: id !< unique identifier matching with the real counterpart
55  integer(I4B) :: container_type !< to identify the actual type of this container
56  character(LENCOMPONENTNAME) :: name !< container name (model, exchange, ...) used in the memory path
57  character(LENCONTEXTNAME) :: vmem_ctx !< prefixes virtual memory located on remote processes
58  logical(LGP) :: is_local !< when true, the physical object resides on the same process. However,
59  !! some of its variables can still be remote
60  logical(LGP) :: is_active !< when true, this container is being synchronized
61  integer(I4B) :: orig_rank !< the global rank of the process which holds the physical data for this container,
62  !< for exchanges this is set to be where model1 sits, or if model1 is local,
63  !< where model2 sits
64  type(stlvecint) :: rcv_ranks !< the ranks of processes, other than orig_rank, having this container active
65  !< (only guaranteed to be complete after synchronization)
66 
67  type(listtype) :: virtual_data_list !< a list with all virtual data items for this container
68  type(vdcelementmaptype), dimension(NR_VDC_ELEMENT_MAPS) :: element_maps !< a list with all element maps
69  type(vdcelementluttype), dimension(NR_VDC_ELEMENT_MAPS) :: element_luts !< lookup tables from remote index to local index
70  contains
71  procedure :: vdc_create
72  generic :: map => map_scalar, map_array1d, map_array2d
73  procedure :: prepare_stage => vdc_prepare_stage
74  procedure :: link_items => vdc_link_items
75  procedure :: set_element_map => vdc_set_element_map
76  procedure :: get_vrt_mem_path => vdc_get_vrt_mem_path
77  procedure :: destroy => vdc_destroy
78  procedure :: set_orig_rank => vdc_set_orig_rank
79  procedure :: get_send_items => vdc_get_send_items
80  procedure :: get_recv_items => vdc_get_recv_items
81  procedure :: get_virtual_data => vdc_get_virtual_data
82  procedure :: print_items
83  ! protected
84  procedure :: set
85  ! private
86  procedure, private :: add_to_list
87  procedure, private :: map_scalar
88  procedure, private :: map_array1d
89  procedure, private :: map_array2d
90  procedure, private :: map_internal
91  procedure, private :: vdc_get_virtual_data
92  procedure, private :: get_items_for_stage
94 
95 contains
96 
97  subroutine vdc_create(this, name, id, is_local)
98  class(virtualdatacontainertype) :: this
99  character(len=*) :: name
100  integer(I4B) :: id
101  logical(LGP) :: is_local
102  ! local
103  integer(I4B) :: i
104 
105  this%name = name
106  this%id = id
107  this%is_local = is_local
108  this%vmem_ctx = 'undefined'
109  this%orig_rank = 0
110  this%is_active = .true.
111  this%container_type = vdc_unknown_type
112 
113  do i = 1, size(this%element_maps)
114  this%element_maps(i)%nr_virt_elems = 0
115  this%element_maps(i)%remote_elem_shift => null()
116  end do
117  do i = 1, size(this%element_luts)
118  this%element_luts(i)%max_remote_idx = 0
119  this%element_luts(i)%remote_to_virtual => null()
120  end do
121 
122  call this%rcv_ranks%init()
123 
124  end subroutine vdc_create
125 
126  !> @brief Init virtual data item, without allocation,
127  !< and store it in this container.
128  subroutine set(this, field, var_name, subcmp_name, map_id, is_local)
129  class(virtualdatacontainertype) :: this
130  class(virtualdatatype), pointer :: field
131  character(len=*) :: var_name
132  character(len=*) :: subcmp_name
133  integer(I4B) :: map_id
134  logical(LGP), optional :: is_local
135 
136  field%is_remote = .not. this%is_local
137  field%map_type = map_id
138  if (present(is_local)) field%is_remote = .not. is_local
139  field%var_name = var_name
140  field%subcmp_name = subcmp_name
141  if (subcmp_name == '') then
142  field%mem_path = create_mem_path(this%name)
143  else
144  field%mem_path = create_mem_path(this%name, subcmp_name)
145  end if
146  field%is_reduced = (field%is_remote .and. field%map_type > 0)
147  field%remote_elem_shift => null()
148  field%remote_to_virtual => null()
149  field%virtual_mt => null()
150  call this%add_to_list(field)
151 
152  end subroutine set
153 
154  subroutine add_to_list(this, virtual_data)
155  class(virtualdatacontainertype) :: this
156  class(virtualdatatype), pointer :: virtual_data
157  ! local
158  class(*), pointer :: vdata_ptr
159 
160  vdata_ptr => virtual_data
161  call this%virtual_data_list%Add(vdata_ptr)
162 
163  end subroutine add_to_list
164 
165  subroutine vdc_prepare_stage(this, stage)
166  use simmodule, only: ustop
167  class(virtualdatacontainertype) :: this
168  integer(I4B) :: stage
169 
170  write (*, *) 'Error: prepare_stage should be overridden'
171  call ustop()
172 
173  end subroutine vdc_prepare_stage
174 
175  !> @brief Link all local data items to memory
176  !<
177  subroutine vdc_link_items(this, stage)
178  class(virtualdatacontainertype) :: this
179  integer(I4B) :: stage
180  ! local
181  integer(I4B) :: i
182  class(*), pointer :: vdi
183 
184  do i = 1, this%virtual_data_list%Count()
185  vdi => this%virtual_data_list%GetItem(i)
186  select type (vdi)
187  class is (virtualdatatype)
188  if (vdi%is_remote) cycle
189  if (vdi%check_stage(stage)) call vdi%link()
190  end select
191  end do
192 
193  end subroutine vdc_link_items
194 
195  !> @brief Add the source indexes associated with map_id
196  !! as a element map to this container, such that
197  !< src_indexes(1:n) = (i_orig_1 - 1, ..., i_orig_n - 1)
198  subroutine vdc_set_element_map(this, src_indexes, map_id)
199  class(virtualdatacontainertype) :: this
200  integer(I4B), dimension(:), pointer, contiguous :: src_indexes
201  integer(I4B) :: map_id
202  ! local
203  integer(I4B) :: i, idx_remote, max_remote_idx
204 
205  if (this%element_maps(map_id)%nr_virt_elems > 0) then
206  write (*, *) "Error, VDC element map already set"
207  call ustop()
208  end if
209 
210  this%element_maps(map_id)%nr_virt_elems = size(src_indexes)
211  allocate (this%element_maps(map_id)%remote_elem_shift(size(src_indexes)))
212  do i = 1, size(src_indexes)
213  this%element_maps(map_id)%remote_elem_shift(i) = src_indexes(i) - 1
214  end do
215 
216  max_remote_idx = maxval(src_indexes)
217  this%element_luts(map_id)%max_remote_idx = max_remote_idx
218  allocate (this%element_luts(map_id)%remote_to_virtual(max_remote_idx))
219  do i = 1, max_remote_idx
220  this%element_luts(map_id)%remote_to_virtual(i) = -1
221  end do
222  do i = 1, size(src_indexes)
223  idx_remote = src_indexes(i)
224  this%element_luts(map_id)%remote_to_virtual(idx_remote) = i
225  end do
226 
227  end subroutine vdc_set_element_map
228 
229  subroutine map_scalar(this, vd, stages)
230  class(virtualdatacontainertype) :: this
231  class(virtualdatatype), pointer :: vd
232  integer(I4B), dimension(:) :: stages
233 
234  call this%map_internal(vd, (/0/), stages)
235 
236  end subroutine map_scalar
237 
238  subroutine map_array1d(this, vd, nrow, stages)
239  class(virtualdatacontainertype) :: this
240  class(virtualdatatype), pointer :: vd
241  integer(I4B) :: nrow
242  integer(I4B), dimension(:) :: stages
243 
244  call this%map_internal(vd, (/nrow/), stages)
245 
246  end subroutine map_array1d
247 
248  subroutine map_array2d(this, vd, ncol, nrow, stages)
249  class(virtualdatacontainertype) :: this
250  class(virtualdatatype), pointer :: vd
251  integer(I4B) :: ncol
252  integer(I4B) :: nrow
253  integer(I4B), dimension(:) :: stages
254 
255  call this%map_internal(vd, (/ncol, nrow/), stages)
256 
257  end subroutine map_array2d
258 
259  subroutine map_internal(this, vd, shape, stages)
260  class(virtualdatacontainertype) :: this
261  class(virtualdatatype), pointer :: vd
262  integer(I4B), dimension(:) :: shape
263  integer(I4B), dimension(:) :: stages
264  ! local
265  character(len=LENMEMPATH) :: vm_pth
266  logical(LGP) :: found
267 
268  vd%sync_stages = stages
269  if (vd%is_remote) then
270  ! create new virtual memory item
271  vm_pth = this%get_vrt_mem_path(vd%var_name, vd%subcmp_name)
272  call vd%vm_allocate(vd%var_name, vm_pth, shape)
273  call get_from_memorystore(vd%var_name, vm_pth, vd%virtual_mt, found)
274  if (vd%map_type > 0) then
275  vd%remote_to_virtual => this%element_luts(vd%map_type)%remote_to_virtual
276  vd%remote_elem_shift => this%element_maps(vd%map_type)%remote_elem_shift
277  end if
278  end if
279 
280  end subroutine map_internal
281 
282  !> @brief Get indexes of virtual data items to be
283  !< sent for a given stage and rank
284  subroutine vdc_get_send_items(this, stg, rank, vi)
285  class(virtualdatacontainertype) :: this
286  integer(I4B) :: stg
287  integer(I4B) :: rank
288  type(stlvecint) :: vi
289 
290  call this%get_items_for_stage(stg, vi)
291 
292  end subroutine vdc_get_send_items
293 
294  !> @brief Get indexes of virtual data items to be
295  !< received for a given stage and rank
296  subroutine vdc_get_recv_items(this, stg, rank, vi)
297  class(virtualdatacontainertype) :: this
298  integer(I4B) :: stg
299  integer(I4B) :: rank
300  type(stlvecint) :: vi
301 
302  call this%get_items_for_stage(stg, vi)
303 
304  end subroutine vdc_get_recv_items
305 
306  subroutine get_items_for_stage(this, stage, virtual_items)
307  class(virtualdatacontainertype) :: this
308  integer(I4B) :: stage
309  type(stlvecint) :: virtual_items
310  ! local
311  integer(I4B) :: i
312  class(*), pointer :: obj_ptr
313 
314  do i = 1, this%virtual_data_list%Count()
315  obj_ptr => this%virtual_data_list%GetItem(i)
316  select type (obj_ptr)
317  class is (virtualdatatype)
318  if (.not. obj_ptr%check_stage(stage)) cycle
319  call virtual_items%push_back(i)
320  end select
321  end do
322 
323  end subroutine get_items_for_stage
324 
325  subroutine print_items(this, imon, items)
326  class(virtualdatacontainertype) :: this
327  integer(I4B) :: imon
328  type(stlvecint) :: items
329  ! local
330  integer(I4B) :: i
331  class(virtualdatatype), pointer :: vdi
332 
333  write (imon, *) "=====> items"
334  do i = 1, items%size
335  vdi => get_virtual_data_from_list(this%virtual_data_list, items%at(i))
336  write (imon, *) vdi%var_name, ":", vdi%mem_path
337  end do
338  if (items%size == 0) then
339  write (imon, *) "... empty ...", this%name
340  end if
341  write (imon, *) "<===== items"
342 
343  end subroutine print_items
344 
345  !> @brief Get virtual memory path for a certain variable
346  !<
347  function vdc_get_vrt_mem_path(this, var_name, subcomp_name) result(vrt_path)
348  class(virtualdatacontainertype) :: this
349  character(len=*) :: var_name
350  character(len=*) :: subcomp_name
351  character(len=LENMEMPATH) :: vrt_path
352  ! local
353  class(virtualdatatype), pointer :: vdi
354 
355  vdi => this%vdc_get_virtual_data(var_name, subcomp_name)
356  if (vdi%is_remote) then
357  if (subcomp_name == '') then
358  vrt_path = create_mem_path(this%name, context=this%vmem_ctx)
359  else
360  vrt_path = create_mem_path(this%name, subcomp_name, context=this%vmem_ctx)
361  end if
362  else
363  if (subcomp_name == '') then
364  vrt_path = create_mem_path(this%name)
365  else
366  vrt_path = create_mem_path(this%name, subcomp_name)
367  end if
368  end if
369 
370  end function vdc_get_vrt_mem_path
371 
372  function vdc_get_virtual_data(this, var_name, subcomp_name) result(virtual_data)
373  use simmodule, only: ustop
374  class(virtualdatacontainertype) :: this
375  character(len=*) :: var_name
376  character(len=*) :: subcomp_name
377  class(virtualdatatype), pointer :: virtual_data
378  ! local
379  integer(I4B) :: i
380  class(virtualdatatype), pointer :: vd
381 
382  virtual_data => null()
383  do i = 1, this%virtual_data_list%Count()
384  vd => get_virtual_data_from_list(this%virtual_data_list, i)
385  if (vd%var_name == var_name .and. &
386  vd%subcmp_name == subcomp_name) then
387  virtual_data => vd
388  return
389  end if
390  end do
391 
392  write (*, *) 'Error: unknown virtual variable ', var_name, ' ', subcomp_name
393  call ustop()
394 
395  end function vdc_get_virtual_data
396 
397  subroutine vdc_destroy(this)
398  class(virtualdatacontainertype) :: this
399  ! local
400  integer(I4B) :: i
401  class(*), pointer :: obj
402 
403  call this%rcv_ranks%destroy()
404 
405  do i = 1, size(this%element_maps)
406  if (associated(this%element_maps(i)%remote_elem_shift)) then
407  deallocate (this%element_maps(i)%remote_elem_shift)
408  end if
409  end do
410  do i = 1, size(this%element_luts)
411  if (associated(this%element_luts(i)%remote_to_virtual)) then
412  deallocate (this%element_luts(i)%remote_to_virtual)
413  end if
414  end do
415 
416  do i = 1, this%virtual_data_list%Count()
417  obj => this%virtual_data_list%GetItem(i)
418  select type (obj)
419  class is (virtualdatatype)
420  if (associated(obj%virtual_mt)) then
421  call obj%vm_deallocate()
422  end if
423  end select
424  end do
425  call this%virtual_data_list%Clear()
426 
427  end subroutine vdc_destroy
428 
429  subroutine vdc_set_orig_rank(this, rank)
430  class(virtualdatacontainertype) :: this
431  integer(I4B) :: rank
432 
433  this%orig_rank = rank
434  write (this%vmem_ctx, '(a,i0,a)') '__P', rank, '__'
435 
436  end subroutine vdc_set_orig_rank
437 
438  function get_vdc_from_list(list, idx) result(vdc)
439  type(listtype) :: list
440  integer(I4B) :: idx
441  class(virtualdatacontainertype), pointer :: vdc
442  ! local
443  class(*), pointer :: obj_ptr
444 
445  vdc => null()
446  obj_ptr => list%GetItem(idx)
447  select type (obj_ptr)
448  class is (virtualdatacontainertype)
449  vdc => obj_ptr
450  end select
451 
452  end function get_vdc_from_list
453 
454  !> @ Converts a virtual container type to its string representation
455  !<
456  function vdc_type_to_str(cntr_type) result(cntr_str)
457  integer(I4B) :: cntr_type
458  character(len=24) :: cntr_str
459 
460  if (cntr_type == vdc_unknown_type) then; cntr_str = "unknown"
461  else if (cntr_type == vdc_gwfmodel_type) then; cntr_str = "GWF Model"
462  else if (cntr_type == vdc_gwtmodel_type) then; cntr_str = "GWT Model"
463  else if (cntr_type == vdc_gwemodel_type) then; cntr_str = "GWE Model"
464  else if (cntr_type == vdc_gwfexg_type) then; cntr_str = "GWF Exchange"
465  else if (cntr_type == vdc_gwtexg_type) then; cntr_str = "GWT Exchange"
466  else if (cntr_type == vdc_gweexg_type) then; cntr_str = "GWE Exchange"
467  else if (cntr_type == vdc_gwfmvr_type) then; cntr_str = "GWF Mover"
468  else if (cntr_type == vdc_gwtmvt_type) then; cntr_str = "GWT Mover"
469  else if (cntr_type == vdc_gwemve_type) then; cntr_str = "GWE Mover"
470  else; cntr_str = "Undefined"
471  end if
472 
473  end function vdc_type_to_str
474 
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lencomponentname
maximum length of a component name
Definition: Constants.f90:18
integer(i4b), parameter lencontextname
maximum length of a memory manager context
Definition: Constants.f90:19
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
subroutine, public get_from_memorystore(name, mem_path, mt, found, check)
@ brief Get a memory type entry from the memory list
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public ustop(stopmess, ioutlocal)
Stop the simulation.
Definition: Sim.f90:312
subroutine destroy(this)
Definition: STLVecInt.f90:183
class(virtualdatatype) function, pointer, public get_virtual_data_from_list(list, idx)
class(virtualdatatype) function, pointer vdc_get_virtual_data(this, var_name, subcomp_name)
integer(i4b), parameter, public vdc_gwtmodel_type
subroutine vdc_create(this, name, id, is_local)
integer(i4b), parameter, public vdc_gwtmvt_type
subroutine map_internal(this, vd, shape, stages)
subroutine set(this, field, var_name, subcmp_name, map_id, is_local)
Init virtual data item, without allocation,.
integer(i4b), parameter, public vdc_gwemodel_type
character(len=24) function, public vdc_type_to_str(cntr_type)
@ Converts a virtual container type to its string representation
subroutine map_array2d(this, vd, ncol, nrow, stages)
subroutine vdc_link_items(this, stage)
Link all local data items to memory.
subroutine vdc_get_recv_items(this, stg, rank, vi)
Get indexes of virtual data items to be.
subroutine map_scalar(this, vd, stages)
subroutine vdc_set_element_map(this, src_indexes, map_id)
Add the source indexes associated with map_id as a element map to this container, such that.
subroutine add_to_list(this, virtual_data)
integer(i4b), parameter, public vdc_gwfmvr_type
integer(i4b), parameter, public vdc_unknown_type
character(len=lenmempath) function vdc_get_vrt_mem_path(this, var_name, subcomp_name)
Get virtual memory path for a certain variable.
subroutine vdc_set_orig_rank(this, rank)
integer(i4b), parameter, public vdc_gwemve_type
of VdcPtrType instead.
subroutine vdc_prepare_stage(this, stage)
subroutine vdc_get_send_items(this, stg, rank, vi)
Get indexes of virtual data items to be.
integer(i4b), parameter, public vdc_gwfmodel_type
subroutine print_items(this, imon, items)
integer(i4b), parameter, public vdc_gwtexg_type
class(virtualdatacontainertype) function, pointer, public get_vdc_from_list(list, idx)
subroutine map_array1d(this, vd, nrow, stages)
integer(i4b), parameter, public vdc_gwfexg_type
integer(i4b), parameter, public vdc_gweexg_type
subroutine get_items_for_stage(this, stage, virtual_items)
A generic heterogeneous doubly-linked list.
Definition: List.f90:14
This is a generic data structure to virtualize pieces of memory in 2 distinct ways:
Definition: VirtualBase.f90:35
Wrapper for virtual data containers.