MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
meshdisvmodelmodule Module Reference

This module contains the MeshDisvModelModule. More...

Data Types

type  mesh2ddisvexporttype
 

Functions/Subroutines

subroutine disv_export_init (this, modelname, modeltype, modelfname, nc_fname, disenum, nctype, iout)
 netcdf export disv init More...
 
subroutine disv_export_destroy (this)
 netcdf export disv destroy More...
 
subroutine df (this)
 netcdf export define More...
 
subroutine step (this)
 netcdf export step More...
 
subroutine package_step (this, export_pkg)
 netcdf export package dynamic input More...
 
subroutine export_input_array (this, pkgtype, pkgname, mempath, idt)
 netcdf export an input array More...
 
subroutine define_dim (this)
 netcdf export define dimensions More...
 
subroutine add_mesh_data (this)
 netcdf export add mesh information More...
 
subroutine nc_export_int1d (p_mem, ncid, dim_ids, var_ids, disv, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, iper, nc_fname)
 netcdf export 1D integer array More...
 
subroutine nc_export_int2d (p_mem, ncid, dim_ids, var_ids, disv, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 2D integer array More...
 
subroutine nc_export_dbl1d (p_mem, ncid, dim_ids, var_ids, disv, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, iper, iaux, nc_fname)
 netcdf export 1D double array More...
 
subroutine nc_export_dbl2d (p_mem, ncid, dim_ids, var_ids, disv, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 2D double array More...
 

Detailed Description

This module defines UGRID layered mesh compliant netcdf export type for DISV models. It is dependent on netcdf libraries.

Function/Subroutine Documentation

◆ add_mesh_data()

subroutine meshdisvmodelmodule::add_mesh_data ( class(mesh2ddisvexporttype), intent(inout)  this)
private

Definition at line 400 of file DisvNCMesh.f90.

402  class(Mesh2dDisvExportType), intent(inout) :: this
403  integer(I4B), dimension(:), contiguous, pointer :: icell2d => null()
404  integer(I4B), dimension(:), contiguous, pointer :: ncvert => null()
405  integer(I4B), dimension(:), contiguous, pointer :: icvert => null()
406  real(DP), dimension(:), contiguous, pointer :: cell_x => null()
407  real(DP), dimension(:), contiguous, pointer :: cell_y => null()
408  real(DP), dimension(:), contiguous, pointer :: vert_x => null()
409  real(DP), dimension(:), contiguous, pointer :: vert_y => null()
410  real(DP), dimension(:), contiguous, pointer :: cell_xt => null()
411  real(DP), dimension(:), contiguous, pointer :: cell_yt => null()
412  real(DP), dimension(:), contiguous, pointer :: vert_xt => null()
413  real(DP), dimension(:), contiguous, pointer :: vert_yt => null()
414  real(DP) :: x_transform, y_transform
415  integer(I4B) :: n, m, idx, cnt, iv, maxvert
416  integer(I4B), dimension(:), allocatable :: verts
417  real(DP), dimension(:), allocatable :: bnds
418  integer(I4B) :: istop
419 
420  ! set pointers to input context
421  call mem_setptr(icell2d, 'ICELL2D', this%dis_mempath)
422  call mem_setptr(ncvert, 'NCVERT', this%dis_mempath)
423  call mem_setptr(icvert, 'ICVERT', this%dis_mempath)
424  call mem_setptr(cell_x, 'XC', this%dis_mempath)
425  call mem_setptr(cell_y, 'YC', this%dis_mempath)
426  call mem_setptr(vert_x, 'XV', this%dis_mempath)
427  call mem_setptr(vert_y, 'YV', this%dis_mempath)
428 
429  ! allocate x, y transform arrays
430  allocate (cell_xt(size(cell_x)))
431  allocate (cell_yt(size(cell_y)))
432  allocate (vert_xt(size(vert_x)))
433  allocate (vert_yt(size(vert_y)))
434 
435  ! set mesh container variable value to 1
436  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh, 1), &
437  this%nc_fname)
438 
439  ! transform vert x and y
440  do n = 1, size(vert_x)
441  call dis_transform_xy(vert_x(n), vert_y(n), &
442  this%disv%xorigin, &
443  this%disv%yorigin, &
444  this%disv%angrot, &
445  x_transform, y_transform)
446  vert_xt(n) = x_transform
447  vert_yt(n) = y_transform
448  end do
449 
450  ! transform cell x and y
451  do n = 1, size(cell_x)
452  call dis_transform_xy(cell_x(n), cell_y(n), &
453  this%disv%xorigin, &
454  this%disv%yorigin, &
455  this%disv%angrot, &
456  x_transform, y_transform)
457  cell_xt(n) = x_transform
458  cell_yt(n) = y_transform
459  end do
460 
461  ! write node_x and node_y arrays to netcdf file
462  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_x, &
463  vert_xt), this%nc_fname)
464  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_y, &
465  vert_yt), this%nc_fname)
466 
467  ! write face_x and face_y arrays to netcdf file
468  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_x, &
469  cell_xt), this%nc_fname)
470  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_y, &
471  cell_yt), this%nc_fname)
472 
473  ! initialize max vertices required to define cell
474  maxvert = maxval(ncvert)
475 
476  ! allocate temporary arrays
477  allocate (verts(maxvert))
478  allocate (bnds(maxvert))
479 
480  ! set face nodes array
481  cnt = 0
482  do n = 1, size(ncvert)
483  verts = nf90_fill_int
484  idx = cnt + ncvert(n)
485  iv = 0
486  istop = cnt + 1
487  do m = idx, istop, -1
488  cnt = cnt + 1
489  iv = iv + 1
490  verts(iv) = icvert(m)
491  end do
492 
493  ! write face nodes array to netcdf file
494  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_nodes, &
495  verts, start=(/1, n/), &
496  count=(/maxvert, 1/)), &
497  this%nc_fname)
498 
499  ! set face y bounds array
500  bnds = nf90_fill_double
501  do m = 1, size(bnds)
502  if (verts(m) /= nf90_fill_int) then
503  bnds(m) = vert_yt(verts(m))
504  end if
505  ! write face y bounds array to netcdf file
506  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_ybnds, &
507  bnds, start=(/1, n/), &
508  count=(/maxvert, 1/)), &
509  this%nc_fname)
510  end do
511 
512  ! set face x bounds array
513  bnds = nf90_fill_double
514  do m = 1, size(bnds)
515  if (verts(m) /= nf90_fill_int) then
516  bnds(m) = vert_xt(verts(m))
517  end if
518  ! write face x bounds array to netcdf file
519  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_xbnds, &
520  bnds, start=(/1, n/), &
521  count=(/maxvert, 1/)), &
522  this%nc_fname)
523  end do
524  end do
525 
526  ! cleanup
527  deallocate (bnds)
528  deallocate (verts)
529  deallocate (cell_xt)
530  deallocate (cell_yt)
531  deallocate (vert_xt)
532  deallocate (vert_yt)
subroutine, public dis_transform_xy(x, y, xorigin, yorigin, angrot, xglo, yglo)
Get global (x, y) coordinates from cell-local coordinates.
Here is the call graph for this function:

◆ define_dim()

subroutine meshdisvmodelmodule::define_dim ( class(mesh2ddisvexporttype), intent(inout)  this)
private

Definition at line 359 of file DisvNCMesh.f90.

360  class(Mesh2dDisvExportType), intent(inout) :: this
361  integer(I4B), dimension(:), contiguous, pointer :: ncvert
362 
363  ! set pointers to input context
364  call mem_setptr(ncvert, 'NCVERT', this%dis_mempath)
365 
366  ! time
367  call nf_verify(nf90_def_dim(this%ncid, 'time', this%totnstp, &
368  this%dim_ids%time), this%nc_fname)
369  call nf_verify(nf90_def_var(this%ncid, 'time', nf90_double, &
370  this%dim_ids%time, this%var_ids%time), &
371  this%nc_fname)
372  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'calendar', &
373  'standard'), this%nc_fname)
374  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'units', &
375  this%datetime), this%nc_fname)
376  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'axis', 'T'), &
377  this%nc_fname)
378  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'standard_name', &
379  'time'), this%nc_fname)
380  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'long_name', &
381  'time'), this%nc_fname)
382 
383  ! mesh
384  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_node', this%disv%nvert, &
385  this%dim_ids%nmesh_node), this%nc_fname)
386  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_face', this%disv%ncpl, &
387  this%dim_ids%nmesh_face), this%nc_fname)
388  call nf_verify(nf90_def_dim(this%ncid, 'max_nmesh_face_nodes', &
389  maxval(ncvert), &
390  this%dim_ids%max_nmesh_face_nodes), &
391  this%nc_fname)
392 
393  ! ncpl, nlay
394  call nf_verify(nf90_def_dim(this%ncid, 'nlay', this%disv%nlay, &
395  this%dim_ids%nlay), this%nc_fname)
Here is the call graph for this function:

◆ df()

subroutine meshdisvmodelmodule::df ( class(mesh2ddisvexporttype), intent(inout)  this)
private

Definition at line 84 of file DisvNCMesh.f90.

85  use constantsmodule, only: mvalidate
86  use simvariablesmodule, only: isim_mode
87  class(Mesh2dDisvExportType), intent(inout) :: this
88  ! put root group file scope attributes
89  call this%add_global_att()
90  ! define root group dimensions and coordinate variables
91  call this%define_dim()
92  ! define mesh variables
93  call this%create_mesh()
94  if (isim_mode /= mvalidate) then
95  ! define the dependent variable
96  call this%define_dependent()
97  end if
98  ! define period input arrays
99  call this%df_export()
100  ! exit define mode
101  call nf_verify(nf90_enddef(this%ncid), this%nc_fname)
102  ! create mesh
103  call this%add_mesh_data()
104  ! define and set package input griddata
105  call this%add_pkg_data()
106  ! define and set gridmap variable
107  call this%define_gridmap()
108  ! synchronize file
109  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
This module contains simulation constants.
Definition: Constants.f90:9
@ mvalidate
validation mode - do not run time steps
Definition: Constants.f90:205
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) isim_mode
simulation mode
Here is the call graph for this function:

◆ disv_export_destroy()

subroutine meshdisvmodelmodule::disv_export_destroy ( class(mesh2ddisvexporttype), intent(inout)  this)

Definition at line 74 of file DisvNCMesh.f90.

75  class(Mesh2dDisvExportType), intent(inout) :: this
76  deallocate (this%var_ids%dependent)
77  ! destroy base class
78  call this%mesh_destroy()
79  call this%NCModelExportType%destroy()

◆ disv_export_init()

subroutine meshdisvmodelmodule::disv_export_init ( class(mesh2ddisvexporttype), intent(inout)  this,
character(len=*), intent(in)  modelname,
character(len=*), intent(in)  modeltype,
character(len=*), intent(in)  modelfname,
character(len=*), intent(in)  nc_fname,
integer(i4b), intent(in)  disenum,
integer(i4b), intent(in)  nctype,
integer(i4b), intent(in)  iout 
)
private

Definition at line 48 of file DisvNCMesh.f90.

51  class(Mesh2dDisvExportType), intent(inout) :: this
52  character(len=*), intent(in) :: modelname
53  character(len=*), intent(in) :: modeltype
54  character(len=*), intent(in) :: modelfname
55  character(len=*), intent(in) :: nc_fname
56  integer(I4B), intent(in) :: disenum
57  integer(I4B), intent(in) :: nctype
58  integer(I4B), intent(in) :: iout
59 
60  ! set nlay
61  this%nlay = this%disv%nlay
62 
63  ! allocate var_id arrays
64  allocate (this%var_ids%dependent(this%nlay))
65  allocate (this%var_ids%export(this%nlay))
66 
67  ! initialize base class
68  call this%mesh_init(modelname, modeltype, modelfname, nc_fname, disenum, &
69  nctype, this%disv%lenuni, iout)

◆ export_input_array()

subroutine meshdisvmodelmodule::export_input_array ( class(mesh2ddisvexporttype), intent(inout)  this,
character(len=*), intent(in)  pkgtype,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  mempath,
type(inputparamdefinitiontype), intent(in), pointer  idt 
)

Definition at line 308 of file DisvNCMesh.f90.

309  class(Mesh2dDisvExportType), intent(inout) :: this
310  character(len=*), intent(in) :: pkgtype
311  character(len=*), intent(in) :: pkgname
312  character(len=*), intent(in) :: mempath
313  type(InputParamDefinitionType), pointer, intent(in) :: idt
314  integer(I4B), dimension(:), pointer, contiguous :: int1d
315  integer(I4B), dimension(:, :), pointer, contiguous :: int2d
316  real(DP), dimension(:), pointer, contiguous :: dbl1d
317  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
318  character(len=LINELENGTH) :: nc_tag
319  integer(I4B) :: iper, iaux
320 
321  iper = 0
322  iaux = 0
323 
324  ! set variable input tag
325  nc_tag = this%input_attribute(pkgname, idt)
326 
327  select case (idt%datatype)
328  case ('INTEGER1D')
329  call mem_setptr(int1d, idt%mf6varname, mempath)
330  call nc_export_int1d(int1d, this%ncid, this%dim_ids, this%var_ids, &
331  this%disv, idt, mempath, nc_tag, pkgname, &
332  this%gridmap_name, this%deflate, this%shuffle, &
333  this%chunk_face, iper, this%nc_fname)
334  case ('INTEGER2D')
335  call mem_setptr(int2d, idt%mf6varname, mempath)
336  call nc_export_int2d(int2d, this%ncid, this%dim_ids, this%var_ids, &
337  this%disv, idt, mempath, nc_tag, pkgname, &
338  this%gridmap_name, this%deflate, this%shuffle, &
339  this%chunk_face, this%nc_fname)
340  case ('DOUBLE1D')
341  call mem_setptr(dbl1d, idt%mf6varname, mempath)
342  call nc_export_dbl1d(dbl1d, this%ncid, this%dim_ids, this%var_ids, &
343  this%disv, idt, mempath, nc_tag, pkgname, &
344  this%gridmap_name, this%deflate, this%shuffle, &
345  this%chunk_face, iper, iaux, this%nc_fname)
346  case ('DOUBLE2D')
347  call mem_setptr(dbl2d, idt%mf6varname, mempath)
348  call nc_export_dbl2d(dbl2d, this%ncid, this%dim_ids, this%var_ids, &
349  this%disv, idt, mempath, nc_tag, pkgname, &
350  this%gridmap_name, this%deflate, this%shuffle, &
351  this%chunk_face, this%nc_fname)
352  case default
353  ! no-op, no other datatypes exported
354  end select
Here is the call graph for this function:

◆ nc_export_dbl1d()

subroutine meshdisvmodelmodule::nc_export_dbl1d ( real(dp), dimension(:), intent(in), pointer, contiguous  p_mem,
integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(disvtype), intent(in), pointer  disv,
type(inputparamdefinitiontype), pointer  idt,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  gridmap_name,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
integer(i4b), intent(in)  iper,
integer(i4b), intent(in)  iaux,
character(len=*), intent(in)  nc_fname 
)
private

Definition at line 728 of file DisvNCMesh.f90.

731  use netcdfcommonmodule, only: ixstp
732  real(DP), dimension(:), pointer, contiguous, intent(in) :: p_mem
733  integer(I4B), intent(in) :: ncid
734  type(MeshNCDimIdType), intent(inout) :: dim_ids
735  type(MeshNCVarIdType), intent(inout) :: var_ids
736  type(DisvType), pointer, intent(in) :: disv
737  type(InputParamDefinitionType), pointer :: idt
738  character(len=*), intent(in) :: mempath
739  character(len=*), intent(in) :: nc_tag
740  character(len=*), intent(in) :: pkgname
741  character(len=*), intent(in) :: gridmap_name
742  integer(I4B), intent(in) :: deflate
743  integer(I4B), intent(in) :: shuffle
744  integer(I4B), intent(in) :: chunk_face
745  integer(I4B), intent(in) :: iper
746  integer(I4B), intent(in) :: iaux
747  character(len=*), intent(in) :: nc_fname
748  real(DP), dimension(:), pointer, contiguous :: dbl1d
749  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
750  integer(I4B) :: axis_sz, k, istp
751  integer(I4B), dimension(:), allocatable :: var_id
752  character(len=LINELENGTH) :: longname, varname
753 
754  if (idt%shape == 'NCPL' .or. &
755  idt%shape == 'NAUX NCPL') then
756 
757  if (iper == 0) then
758  ! set names
759  varname = export_varname(pkgname, idt%tagname, mempath, &
760  iaux=iaux)
761  longname = export_longname(idt%longname, pkgname, idt%tagname, &
762  mempath, iaux=iaux)
763 
764  allocate (var_id(1))
765  axis_sz = dim_ids%nmesh_face
766 
767  ! reenter define mode and create variable
768  call nf_verify(nf90_redef(ncid), nc_fname)
769  call nf_verify(nf90_def_var(ncid, varname, nf90_double, &
770  (/axis_sz/), var_id(1)), &
771  nc_fname)
772 
773  ! apply chunking parameters
774  call ncvar_chunk(ncid, var_id(1), chunk_face, nc_fname)
775  ! deflate and shuffle
776  call ncvar_deflate(ncid, var_id(1), deflate, shuffle, nc_fname)
777 
778  ! put attr
779  call nf_verify(nf90_put_att(ncid, var_id(1), '_FillValue', &
780  (/nf90_fill_double/)), nc_fname)
781  call nf_verify(nf90_put_att(ncid, var_id(1), 'long_name', &
782  longname), nc_fname)
783 
784  ! add grid mapping and mf6 attr
785  call ncvar_gridmap(ncid, var_id(1), gridmap_name, nc_fname)
786  call ncvar_mf6attr(ncid, var_id(1), 0, iaux, nc_tag, nc_fname)
787 
788  ! exit define mode and write data
789  call nf_verify(nf90_enddef(ncid), nc_fname)
790  call nf_verify(nf90_put_var(ncid, var_id(1), p_mem), &
791  nc_fname)
792  else
793  ! timeseries
794  istp = ixstp()
795  call nf_verify(nf90_put_var(ncid, &
796  var_ids%export(1), p_mem, &
797  start=(/1, istp/), &
798  count=(/disv%ncpl, 1/)), nc_fname)
799  end if
800 
801  else
802 
803  dbl2d(1:disv%ncpl, 1:disv%nlay) => p_mem(1:disv%nodesuser)
804 
805  if (iper == 0) then
806  allocate (var_id(disv%nlay))
807 
808  ! reenter define mode and create variable
809  call nf_verify(nf90_redef(ncid), nc_fname)
810  do k = 1, disv%nlay
811  ! set names
812  varname = export_varname(pkgname, idt%tagname, mempath, layer=k, &
813  iaux=iaux)
814  longname = export_longname(idt%longname, pkgname, idt%tagname, &
815  mempath, layer=k, iaux=iaux)
816 
817  call nf_verify(nf90_def_var(ncid, varname, nf90_double, &
818  (/dim_ids%nmesh_face/), var_id(k)), &
819  nc_fname)
820 
821  ! apply chunking parameters
822  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
823  ! deflate and shuffle
824  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
825 
826  ! put attr
827  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
828  (/nf90_fill_double/)), nc_fname)
829  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
830  longname), nc_fname)
831 
832  ! add grid mapping and mf6 attr
833  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
834  call ncvar_mf6attr(ncid, var_id(k), k, iaux, nc_tag, nc_fname)
835  end do
836 
837  ! exit define mode and write data
838  call nf_verify(nf90_enddef(ncid), nc_fname)
839  do k = 1, disv%nlay
840  call nf_verify(nf90_put_var(ncid, var_id(k), dbl2d(:, k)), nc_fname)
841  end do
842 
843  ! cleanup
844  deallocate (var_id)
845  else
846  ! timeseries, add period data
847  istp = ixstp()
848  do k = 1, disv%nlay
849  dbl1d(1:disv%ncpl) => dbl2d(:, k)
850  call nf_verify(nf90_put_var(ncid, &
851  var_ids%export(k), dbl1d, &
852  start=(/1, istp/), &
853  count=(/disv%ncpl, 1/)), nc_fname)
854  end do
855  end if
856  end if
This module contains the NetCDFCommonModule.
Definition: NetCDFCommon.f90:6
integer(i4b) function, public ixstp()
step index for timeseries data
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_dbl2d()

subroutine meshdisvmodelmodule::nc_export_dbl2d ( real(dp), dimension(:, :), intent(in), pointer, contiguous  p_mem,
integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(disvtype), intent(in), pointer  disv,
type(inputparamdefinitiontype), pointer  idt,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  gridmap_name,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)

Definition at line 861 of file DisvNCMesh.f90.

864  real(DP), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
865  integer(I4B), intent(in) :: ncid
866  type(MeshNCDimIdType), intent(inout) :: dim_ids
867  type(MeshNCVarIdType), intent(inout) :: var_ids
868  type(DisvType), pointer, intent(in) :: disv
869  type(InputParamDefinitionType), pointer :: idt
870  character(len=*), intent(in) :: mempath
871  character(len=*), intent(in) :: nc_tag
872  character(len=*), intent(in) :: pkgname
873  character(len=*), intent(in) :: gridmap_name
874  integer(I4B), intent(in) :: deflate
875  integer(I4B), intent(in) :: shuffle
876  integer(I4B), intent(in) :: chunk_face
877  character(len=*), intent(in) :: nc_fname
878  integer(I4B), dimension(:), allocatable :: var_id
879  character(len=LINELENGTH) :: longname, varname
880  integer(I4B) :: k
881 
882  allocate (var_id(disv%nlay))
883 
884  ! reenter define mode and create variable
885  call nf_verify(nf90_redef(ncid), nc_fname)
886  do k = 1, disv%nlay
887  ! set names
888  varname = export_varname(pkgname, idt%tagname, mempath, layer=k)
889  longname = export_longname(idt%longname, pkgname, idt%tagname, &
890  mempath, layer=k)
891 
892  call nf_verify(nf90_def_var(ncid, varname, nf90_double, &
893  (/dim_ids%nmesh_face/), var_id(k)), &
894  nc_fname)
895 
896  ! apply chunking parameters
897  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
898  ! deflate and shuffle
899  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
900 
901  ! put attr
902  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
903  (/nf90_fill_double/)), nc_fname)
904  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
905  longname), nc_fname)
906 
907  ! add grid mapping and mf6 attr
908  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
909  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
910  end do
911 
912  ! exit define mode and write data
913  call nf_verify(nf90_enddef(ncid), nc_fname)
914  do k = 1, disv%nlay
915  call nf_verify(nf90_put_var(ncid, var_id(k), p_mem(:, k)), nc_fname)
916  end do
917 
918  deallocate (var_id)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_int1d()

subroutine meshdisvmodelmodule::nc_export_int1d ( integer(i4b), dimension(:), intent(in), pointer, contiguous  p_mem,
integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(disvtype), intent(in), pointer  disv,
type(inputparamdefinitiontype), pointer  idt,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  gridmap_name,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
integer(i4b), intent(in)  iper,
character(len=*), intent(in)  nc_fname 
)

Definition at line 537 of file DisvNCMesh.f90.

540  use netcdfcommonmodule, only: ixstp
541  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: p_mem
542  integer(I4B), intent(in) :: ncid
543  type(MeshNCDimIdType), intent(inout) :: dim_ids
544  type(MeshNCVarIdType), intent(inout) :: var_ids
545  type(DisvType), pointer, intent(in) :: disv
546  type(InputParamDefinitionType), pointer :: idt
547  character(len=*), intent(in) :: mempath
548  character(len=*), intent(in) :: nc_tag
549  character(len=*), intent(in) :: pkgname
550  character(len=*), intent(in) :: gridmap_name
551  integer(I4B), intent(in) :: deflate
552  integer(I4B), intent(in) :: shuffle
553  integer(I4B), intent(in) :: chunk_face
554  integer(I4B), intent(in) :: iper
555  character(len=*), intent(in) :: nc_fname
556  integer(I4B), dimension(:), pointer, contiguous :: int1d
557  integer(I4B), dimension(:, :), pointer, contiguous :: int2d
558  integer(I4B) :: axis_sz, k, istp
559  integer(I4B), dimension(:), allocatable :: var_id
560  character(len=LINELENGTH) :: longname, varname
561 
562  if (idt%shape == 'NCPL' .or. &
563  idt%shape == 'NAUX NCPL') then
564 
565  if (iper == 0) then
566  ! set names
567  varname = export_varname(pkgname, idt%tagname, mempath)
568  longname = export_longname(idt%longname, pkgname, idt%tagname, mempath)
569 
570  allocate (var_id(1))
571  axis_sz = dim_ids%nmesh_face
572 
573  ! reenter define mode and create variable
574  call nf_verify(nf90_redef(ncid), nc_fname)
575  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
576  (/axis_sz/), var_id(1)), &
577  nc_fname)
578 
579  ! apply chunking parameters
580  call ncvar_chunk(ncid, var_id(1), chunk_face, nc_fname)
581  ! deflate and shuffle
582  call ncvar_deflate(ncid, var_id(1), deflate, shuffle, nc_fname)
583 
584  ! put attr
585  call nf_verify(nf90_put_att(ncid, var_id(1), '_FillValue', &
586  (/nf90_fill_int/)), nc_fname)
587  call nf_verify(nf90_put_att(ncid, var_id(1), 'long_name', &
588  longname), nc_fname)
589 
590  ! add grid mapping and mf6 attr
591  call ncvar_gridmap(ncid, var_id(1), gridmap_name, nc_fname)
592  call ncvar_mf6attr(ncid, var_id(1), 0, 0, nc_tag, nc_fname)
593 
594  ! exit define mode and write data
595  call nf_verify(nf90_enddef(ncid), nc_fname)
596  call nf_verify(nf90_put_var(ncid, var_id(1), p_mem), &
597  nc_fname)
598  else
599  ! timeseries
600  istp = ixstp()
601  call nf_verify(nf90_put_var(ncid, &
602  var_ids%export(1), p_mem, &
603  start=(/1, istp/), &
604  count=(/disv%ncpl, 1/)), nc_fname)
605  end if
606 
607  else
608 
609  int2d(1:disv%ncpl, 1:disv%nlay) => p_mem(1:disv%nodesuser)
610 
611  if (iper == 0) then
612  allocate (var_id(disv%nlay))
613 
614  ! reenter define mode and create variable
615  call nf_verify(nf90_redef(ncid), nc_fname)
616  do k = 1, disv%nlay
617  ! set names
618  varname = export_varname(pkgname, idt%tagname, mempath, layer=k)
619  longname = export_longname(idt%longname, pkgname, idt%tagname, &
620  mempath, layer=k)
621 
622  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
623  (/dim_ids%nmesh_face/), var_id(k)), &
624  nc_fname)
625 
626  ! apply chunking parameters
627  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
628  ! deflate and shuffle
629  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
630 
631  ! put attr
632  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
633  (/nf90_fill_int/)), nc_fname)
634  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
635  longname), nc_fname)
636 
637  ! add grid mapping and mf6 attr
638  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
639  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
640  end do
641 
642  ! exit define mode and write data
643  call nf_verify(nf90_enddef(ncid), nc_fname)
644  do k = 1, disv%nlay
645  call nf_verify(nf90_put_var(ncid, var_id(k), int2d(:, k)), nc_fname)
646  end do
647 
648  ! cleanup
649  deallocate (var_id)
650  else
651  ! timeseries, add period data
652  istp = ixstp()
653  do k = 1, disv%nlay
654  int1d(1:disv%ncpl) => int2d(:, k)
655  call nf_verify(nf90_put_var(ncid, &
656  var_ids%export(k), int1d, &
657  start=(/1, istp/), &
658  count=(/disv%ncpl, 1/)), nc_fname)
659  end do
660  end if
661  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_int2d()

subroutine meshdisvmodelmodule::nc_export_int2d ( integer(i4b), dimension(:, :), intent(in), pointer, contiguous  p_mem,
integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(disvtype), intent(in), pointer  disv,
type(inputparamdefinitiontype), pointer  idt,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  gridmap_name,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)

Definition at line 666 of file DisvNCMesh.f90.

669  integer(I4B), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
670  integer(I4B), intent(in) :: ncid
671  type(MeshNCDimIdType), intent(inout) :: dim_ids
672  type(MeshNCVarIdType), intent(inout) :: var_ids
673  type(DisvType), pointer, intent(in) :: disv
674  type(InputParamDefinitionType), pointer :: idt
675  character(len=*), intent(in) :: mempath
676  character(len=*), intent(in) :: nc_tag
677  character(len=*), intent(in) :: pkgname
678  character(len=*), intent(in) :: gridmap_name
679  integer(I4B), intent(in) :: deflate
680  integer(I4B), intent(in) :: shuffle
681  integer(I4B), intent(in) :: chunk_face
682  character(len=*), intent(in) :: nc_fname
683  integer(I4B), dimension(:), allocatable :: var_id
684  character(len=LINELENGTH) :: longname, varname
685  integer(I4B) :: k
686 
687  allocate (var_id(disv%nlay))
688 
689  ! reenter define mode and create variable
690  call nf_verify(nf90_redef(ncid), nc_fname)
691  do k = 1, disv%nlay
692  ! set names
693  varname = export_varname(pkgname, idt%tagname, mempath, layer=k)
694  longname = export_longname(idt%longname, pkgname, idt%tagname, &
695  mempath, layer=k)
696 
697  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
698  (/dim_ids%nmesh_face/), var_id(k)), &
699  nc_fname)
700 
701  ! apply chunking parameters
702  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
703  ! deflate and shuffle
704  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
705 
706  ! put attr
707  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
708  (/nf90_fill_int/)), nc_fname)
709  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
710  longname), nc_fname)
711 
712  ! add grid mapping and mf6 attr
713  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
714  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
715  end do
716 
717  ! exit define mode and write data
718  call nf_verify(nf90_enddef(ncid), nc_fname)
719  do k = 1, disv%nlay
720  call nf_verify(nf90_put_var(ncid, var_id(k), p_mem(:, k)), nc_fname)
721  end do
722 
723  deallocate (var_id)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ package_step()

subroutine meshdisvmodelmodule::package_step ( class(mesh2ddisvexporttype), intent(inout)  this,
class(exportpackagetype), intent(in), pointer  export_pkg 
)

Definition at line 182 of file DisvNCMesh.f90.

183  use tdismodule, only: kper
186  class(Mesh2dDisvExportType), intent(inout) :: this
187  class(ExportPackageType), pointer, intent(in) :: export_pkg
188  type(InputParamDefinitionType), pointer :: idt
189  integer(I4B), dimension(:), pointer, contiguous :: int1d
190  real(DP), dimension(:), pointer, contiguous :: dbl1d, nodes
191  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
192  character(len=LINELENGTH) :: nc_tag
193  integer(I4B) :: iaux, iparam, nvals
194  integer(I4B) :: k, n
195  integer(I4B), pointer :: nbound
196 
197  ! initialize
198  iaux = 0
199 
200  ! export defined period input
201  do iparam = 1, export_pkg%nparam
202  ! check if variable was read this period
203  if (export_pkg%param_reads(iparam)%invar < 1) cycle
204 
205  ! set input definition
206  idt => &
207  get_param_definition_type(export_pkg%mf6_input%param_dfns, &
208  export_pkg%mf6_input%component_type, &
209  export_pkg%mf6_input%subcomponent_type, &
210  'PERIOD', export_pkg%param_names(iparam), '')
211 
212  ! set variable input tag
213  nc_tag = this%input_attribute(export_pkg%mf6_input%subcomponent_name, &
214  idt)
215 
216  ! export arrays
217  select case (idt%datatype)
218  case ('INTEGER1D')
219  call mem_setptr(int1d, idt%mf6varname, export_pkg%mf6_input%mempath)
220  this%var_ids%export(1) = export_pkg%varids_param(iparam, 1)
221  call nc_export_int1d(int1d, this%ncid, this%dim_ids, this%var_ids, &
222  this%disv, idt, export_pkg%mf6_input%mempath, &
223  nc_tag, export_pkg%mf6_input%subcomponent_name, &
224  this%gridmap_name, this%deflate, this%shuffle, &
225  this%chunk_face, kper, this%nc_fname)
226  case ('DOUBLE1D')
227  call mem_setptr(dbl1d, idt%mf6varname, export_pkg%mf6_input%mempath)
228  select case (idt%shape)
229  case ('NCPL')
230  this%var_ids%export(1) = export_pkg%varids_param(iparam, 1)
231  call nc_export_dbl1d(dbl1d, this%ncid, this%dim_ids, this%var_ids, &
232  this%disv, idt, export_pkg%mf6_input%mempath, &
233  nc_tag, export_pkg%mf6_input%subcomponent_name, &
234  this%gridmap_name, this%deflate, this%shuffle, &
235  this%chunk_face, kper, iaux, this%nc_fname)
236  case ('NODES')
237  nvals = this%disv%nodesuser
238  allocate (nodes(nvals))
239  nodes = dnodata
240  do k = 1, this%disv%nlay
241  this%var_ids%export(k) = export_pkg%varids_param(iparam, k)
242  end do
243  call mem_setptr(dbl1d, idt%mf6varname, export_pkg%mf6_input%mempath)
244  call mem_setptr(int1d, 'NODEULIST', export_pkg%mf6_input%mempath)
245  call mem_setptr(nbound, 'NBOUND', export_pkg%mf6_input%mempath)
246  do n = 1, nbound
247  nodes(int1d(n)) = dbl1d(n)
248  end do
249  call nc_export_dbl1d(nodes, this%ncid, this%dim_ids, this%var_ids, &
250  this%disv, idt, export_pkg%mf6_input%mempath, &
251  nc_tag, export_pkg%mf6_input%subcomponent_name, &
252  this%gridmap_name, this%deflate, this%shuffle, &
253  this%chunk_face, kper, iaux, this%nc_fname)
254  deallocate (nodes)
255  case default
256  end select
257  case ('DOUBLE2D')
258  call mem_setptr(dbl2d, idt%mf6varname, export_pkg%mf6_input%mempath)
259  select case (idt%shape)
260  case ('NAUX NCPL')
261  nvals = this%disv%ncpl
262  allocate (nodes(nvals))
263  do iaux = 1, size(dbl2d, dim=1) !naux
264  this%var_ids%export(1) = export_pkg%varids_aux(iaux, 1)
265  do n = 1, nvals
266  nodes(n) = dbl2d(iaux, n)
267  end do
268  call nc_export_dbl1d(dbl1d, this%ncid, this%dim_ids, this%var_ids, &
269  this%disv, idt, export_pkg%mf6_input%mempath, &
270  nc_tag, export_pkg%mf6_input%subcomponent_name, &
271  this%gridmap_name, this%deflate, this%shuffle, &
272  this%chunk_face, kper, iaux, this%nc_fname)
273  end do
274  deallocate (nodes)
275  case ('NAUX NODES')
276  nvals = this%disv%nodesuser
277  allocate (nodes(nvals))
278  call mem_setptr(int1d, 'NODEULIST', export_pkg%mf6_input%mempath)
279  call mem_setptr(nbound, 'NBOUND', export_pkg%mf6_input%mempath)
280  do iaux = 1, size(dbl2d, dim=1) ! naux
281  nodes = dnodata
282  do k = 1, this%disv%nlay
283  this%var_ids%export(k) = export_pkg%varids_aux(iaux, k)
284  end do
285  do n = 1, nbound
286  nodes(int1d(n)) = dbl2d(iaux, n)
287  end do
288  call nc_export_dbl1d(nodes, this%ncid, this%dim_ids, this%var_ids, &
289  this%disv, idt, export_pkg%mf6_input%mempath, &
290  nc_tag, export_pkg%mf6_input%subcomponent_name, &
291  this%gridmap_name, this%deflate, this%shuffle, &
292  this%chunk_face, kper, iaux, this%nc_fname)
293  end do
294  deallocate (nodes)
295  case default
296  end select
297  case default
298  ! no-op, no other datatypes exported
299  end select
300  end do
301 
302  ! synchronize file
303  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
This module contains the DefinitionSelectModule.
type(inputparamdefinitiontype) function, pointer, public get_param_definition_type(input_definition_types, component_type, subcomponent_type, blockname, tagname, filename)
Return parameter definition.
This module contains the NCModelExportModule.
Definition: NCModel.f90:8
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ step()

subroutine meshdisvmodelmodule::step ( class(mesh2ddisvexporttype), intent(inout)  this)

Definition at line 114 of file DisvNCMesh.f90.

115  use constantsmodule, only: dhnoflo
116  use tdismodule, only: totim
117  use netcdfcommonmodule, only: ixstp
118  class(Mesh2dDisvExportType), intent(inout) :: this
119  real(DP), dimension(:), pointer, contiguous :: dbl1d
120  integer(I4B) :: n, k, nvals, istp
121  integer(I4B), dimension(2) :: dis_shape
122  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
123 
124  ! initialize
125  nullify (dbl1d)
126  nullify (dbl2d)
127 
128  ! set global step index
129  istp = ixstp()
130 
131  dis_shape(1) = this%disv%ncpl
132  dis_shape(2) = this%disv%nlay
133 
134  nvals = product(dis_shape)
135 
136  ! add data to dependent variable
137  if (size(this%disv%nodeuser) < &
138  size(this%disv%nodereduced)) then
139  ! allocate nodereduced size 1d array
140  allocate (dbl1d(size(this%disv%nodereduced)))
141 
142  ! initialize DHNOFLO for non-active cells
143  dbl1d = dhnoflo
144 
145  ! update active cells
146  do n = 1, size(this%disv%nodereduced)
147  if (this%disv%nodereduced(n) > 0) then
148  dbl1d(n) = this%x(this%disv%nodereduced(n))
149  end if
150  end do
151 
152  dbl2d(1:dis_shape(1), 1:dis_shape(2)) => dbl1d(1:nvals)
153  else
154  dbl2d(1:dis_shape(1), 1:dis_shape(2)) => this%x(1:nvals)
155  end if
156 
157  do k = 1, this%disv%nlay
158  ! extend array with step data
159  call nf_verify(nf90_put_var(this%ncid, &
160  this%var_ids%dependent(k), dbl2d(:, k), &
161  start=(/1, istp/), &
162  count=(/this%disv%ncpl, 1/)), &
163  this%nc_fname)
164  end do
165 
166  ! write to time coordinate variable
167  call nf_verify(nf90_put_var(this%ncid, this%var_ids%time, &
168  totim, start=(/istp/)), &
169  this%nc_fname)
170 
171  ! update file
172  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
173 
174  ! cleanup
175  if (associated(dbl1d)) deallocate (dbl1d)
176  nullify (dbl1d)
177  nullify (dbl2d)
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
Here is the call graph for this function: