MODFLOW 6  version 6.8.0.dev0
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)

Definition at line 407 of file DisvNCMesh.f90.

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

367  use constantsmodule, only: mvalidate
368  use simvariablesmodule, only: isim_mode
369  class(Mesh2dDisvExportType), intent(inout) :: this
370  integer(I4B), dimension(:), contiguous, pointer :: ncvert
371 
372  ! set pointers to input context
373  call mem_setptr(ncvert, 'NCVERT', this%dis_mempath)
374 
375  if (isim_mode /= mvalidate .or. this%pkglist%Count() > 0) then
376  ! time
377  call nf_verify(nf90_def_dim(this%ncid, 'time', this%totnstp, &
378  this%dim_ids%time), this%nc_fname)
379  call nf_verify(nf90_def_var(this%ncid, 'time', nf90_double, &
380  this%dim_ids%time, this%var_ids%time), &
381  this%nc_fname)
382  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'calendar', &
383  'standard'), this%nc_fname)
384  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'units', &
385  this%datetime), this%nc_fname)
386  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'axis', 'T'), &
387  this%nc_fname)
388  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'standard_name', &
389  'time'), this%nc_fname)
390  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'long_name', &
391  'time'), this%nc_fname)
392  end if
393 
394  ! mesh
395  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_node', this%disv%nvert, &
396  this%dim_ids%nmesh_node), this%nc_fname)
397  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_face', this%disv%ncpl, &
398  this%dim_ids%nmesh_face), this%nc_fname)
399  call nf_verify(nf90_def_dim(this%ncid, 'max_nmesh_face_nodes', &
400  maxval(ncvert), &
401  this%dim_ids%max_nmesh_face_nodes), &
402  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:

◆ 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 period input arrays
96  call this%df_export()
97  else
98  ! define the dependent variable
99  call this%define_dependent()
100  end if
101  ! exit define mode
102  call nf_verify(nf90_enddef(this%ncid), this%nc_fname)
103  ! create mesh
104  call this%add_mesh_data()
105  if (isim_mode == mvalidate) then
106  ! define and set package input griddata
107  call this%add_pkg_data()
108  end if
109  ! define and set gridmap variable
110  call this%define_gridmap()
111  ! synchronize file
112  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
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 315 of file DisvNCMesh.f90.

316  class(Mesh2dDisvExportType), intent(inout) :: this
317  character(len=*), intent(in) :: pkgtype
318  character(len=*), intent(in) :: pkgname
319  character(len=*), intent(in) :: mempath
320  type(InputParamDefinitionType), pointer, intent(in) :: idt
321  integer(I4B), dimension(:), pointer, contiguous :: int1d
322  integer(I4B), dimension(:, :), pointer, contiguous :: int2d
323  real(DP), dimension(:), pointer, contiguous :: dbl1d
324  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
325  character(len=LINELENGTH) :: nc_tag
326  integer(I4B) :: iper, iaux
327 
328  iper = 0
329  iaux = 0
330 
331  ! set variable input tag
332  nc_tag = this%input_attribute(pkgname, idt)
333 
334  select case (idt%datatype)
335  case ('INTEGER1D')
336  call mem_setptr(int1d, idt%mf6varname, mempath)
337  call nc_export_int1d(int1d, this%ncid, this%dim_ids, this%var_ids, &
338  this%disv, idt, mempath, nc_tag, pkgname, &
339  this%gridmap_name, this%deflate, this%shuffle, &
340  this%chunk_face, iper, this%nc_fname)
341  case ('INTEGER2D')
342  call mem_setptr(int2d, idt%mf6varname, mempath)
343  call nc_export_int2d(int2d, this%ncid, this%dim_ids, this%var_ids, &
344  this%disv, idt, mempath, nc_tag, pkgname, &
345  this%gridmap_name, this%deflate, this%shuffle, &
346  this%chunk_face, this%nc_fname)
347  case ('DOUBLE1D')
348  call mem_setptr(dbl1d, idt%mf6varname, mempath)
349  call nc_export_dbl1d(dbl1d, this%ncid, this%dim_ids, this%var_ids, &
350  this%disv, idt, mempath, nc_tag, pkgname, &
351  this%gridmap_name, this%deflate, this%shuffle, &
352  this%chunk_face, iper, iaux, this%nc_fname)
353  case ('DOUBLE2D')
354  call mem_setptr(dbl2d, idt%mf6varname, mempath)
355  call nc_export_dbl2d(dbl2d, this%ncid, this%dim_ids, this%var_ids, &
356  this%disv, idt, mempath, nc_tag, pkgname, &
357  this%gridmap_name, this%deflate, this%shuffle, &
358  this%chunk_face, this%nc_fname)
359  case default
360  ! no-op, no other datatypes exported
361  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 736 of file DisvNCMesh.f90.

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

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

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

677  integer(I4B), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
678  integer(I4B), intent(in) :: ncid
679  type(MeshNCDimIdType), intent(inout) :: dim_ids
680  type(MeshNCVarIdType), intent(inout) :: var_ids
681  type(DisvType), pointer, intent(in) :: disv
682  type(InputParamDefinitionType), pointer :: idt
683  character(len=*), intent(in) :: mempath
684  character(len=*), intent(in) :: nc_tag
685  character(len=*), intent(in) :: pkgname
686  character(len=*), intent(in) :: gridmap_name
687  integer(I4B), intent(in) :: deflate
688  integer(I4B), intent(in) :: shuffle
689  integer(I4B), intent(in) :: chunk_face
690  character(len=*), intent(in) :: nc_fname
691  integer(I4B), dimension(:), allocatable :: var_id
692  character(len=LINELENGTH) :: longname, varname
693  integer(I4B) :: k
694 
695  allocate (var_id(disv%nlay))
696 
697  ! reenter define mode and create variable
698  call nf_verify(nf90_redef(ncid), nc_fname)
699  do k = 1, disv%nlay
700  ! set names
701  varname = export_varname(pkgname, idt%tagname, mempath, layer=k)
702  longname = export_longname(idt%longname, pkgname, idt%tagname, &
703  mempath, layer=k)
704 
705  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
706  (/dim_ids%nmesh_face/), var_id(k)), &
707  nc_fname)
708 
709  ! apply chunking parameters
710  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
711  ! deflate and shuffle
712  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
713 
714  ! put attr
715  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
716  (/nf90_fill_int/)), nc_fname)
717  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
718  longname), nc_fname)
719 
720  ! add grid mapping and mf6 attr
721  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
722  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
723  end do
724 
725  ! exit define mode and write data
726  call nf_verify(nf90_enddef(ncid), nc_fname)
727  do k = 1, disv%nlay
728  call nf_verify(nf90_put_var(ncid, var_id(k), p_mem(:, k)), nc_fname)
729  end do
730 
731  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 184 of file DisvNCMesh.f90.

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

◆ step()

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

Definition at line 117 of file DisvNCMesh.f90.

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