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 :: ncvert => null()
411  integer(I4B), dimension(:), contiguous, pointer :: icvert => null()
412  real(DP), dimension(:), contiguous, pointer :: cell_xt => null()
413  real(DP), dimension(:), contiguous, pointer :: cell_yt => null()
414  real(DP), dimension(:), contiguous, pointer :: vert_xt => null()
415  real(DP), dimension(:), contiguous, pointer :: vert_yt => null()
416  real(DP) :: x_transform, y_transform
417  integer(I4B) :: n, m, idx, cnt, iv, maxvert
418  integer(I4B), dimension(:), allocatable :: verts
419  real(DP), dimension(:), allocatable :: bnds
420  integer(I4B) :: istop
421 
422  ! set pointers to input context
423  call mem_setptr(ncvert, 'NCVERT', this%dis_mempath)
424  call mem_setptr(icvert, 'ICVERT', this%dis_mempath)
425 
426  ! allocate x, y transform arrays
427  allocate (cell_xt(this%disv%ncpl))
428  allocate (cell_yt(this%disv%ncpl))
429  allocate (vert_xt(this%disv%nvert))
430  allocate (vert_yt(this%disv%nvert))
431 
432  ! set mesh container variable value to 1
433  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh, 1), &
434  this%nc_fname)
435 
436  ! transform vert x and y
437  do n = 1, this%disv%nvert
438  call dis_transform_xy(this%disv%vertices(1, n), this%disv%vertices(2, n), &
439  this%disv%xorigin, &
440  this%disv%yorigin, &
441  this%disv%angrot, &
442  x_transform, y_transform)
443  vert_xt(n) = x_transform
444  vert_yt(n) = y_transform
445  end do
446 
447  ! transform cell x and y
448  do n = 1, this%disv%ncpl
449  call dis_transform_xy(this%disv%cellxy(1, n), this%disv%cellxy(2, n), &
450  this%disv%xorigin, &
451  this%disv%yorigin, &
452  this%disv%angrot, &
453  x_transform, y_transform)
454  cell_xt(n) = x_transform
455  cell_yt(n) = y_transform
456  end do
457 
458  ! write node_x and node_y arrays to netcdf file
459  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_x, &
460  vert_xt), this%nc_fname)
461  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_y, &
462  vert_yt), this%nc_fname)
463 
464  ! write face_x and face_y arrays to netcdf file
465  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_x, &
466  cell_xt), this%nc_fname)
467  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_y, &
468  cell_yt), this%nc_fname)
469 
470  ! initialize max vertices required to define cell
471  maxvert = maxval(ncvert)
472 
473  ! allocate temporary arrays
474  allocate (verts(maxvert))
475  allocate (bnds(maxvert))
476 
477  ! set face nodes array
478  cnt = 0
479  do n = 1, this%disv%ncpl
480  verts = nf90_fill_int
481  idx = cnt + ncvert(n)
482  iv = 0
483  istop = cnt + 1
484  do m = idx, istop, -1
485  cnt = cnt + 1
486  iv = iv + 1
487  verts(iv) = icvert(m)
488  end do
489 
490  ! don't close the cell
491  if (verts(iv) == verts(1)) verts(iv) = nf90_fill_int
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 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 726 of file DisvNCMesh.f90.

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

860  real(DP), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
861  integer(I4B), intent(in) :: ncid
862  type(MeshNCDimIdType), intent(inout) :: dim_ids
863  type(MeshNCVarIdType), intent(inout) :: var_ids
864  type(DisvType), pointer, intent(in) :: disv
865  type(InputParamDefinitionType), pointer :: idt
866  character(len=*), intent(in) :: mempath
867  character(len=*), intent(in) :: nc_tag
868  character(len=*), intent(in) :: pkgname
869  character(len=*), intent(in) :: gridmap_name
870  integer(I4B), intent(in) :: deflate
871  integer(I4B), intent(in) :: shuffle
872  integer(I4B), intent(in) :: chunk_face
873  character(len=*), intent(in) :: nc_fname
874  integer(I4B), dimension(:), allocatable :: var_id
875  character(len=LINELENGTH) :: longname, varname
876  integer(I4B) :: k
877 
878  allocate (var_id(disv%nlay))
879 
880  ! reenter define mode and create variable
881  call nf_verify(nf90_redef(ncid), nc_fname)
882  do k = 1, disv%nlay
883  ! set names
884  varname = export_varname(pkgname, idt%tagname, mempath, layer=k)
885  longname = export_longname(idt%longname, pkgname, idt%tagname, &
886  mempath, layer=k)
887 
888  call nf_verify(nf90_def_var(ncid, varname, nf90_double, &
889  (/dim_ids%nmesh_face/), var_id(k)), &
890  nc_fname)
891 
892  ! apply chunking parameters
893  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
894  ! deflate and shuffle
895  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
896 
897  ! put attr
898  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
899  (/nf90_fill_double/)), nc_fname)
900  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
901  longname), nc_fname)
902 
903  ! add grid mapping and mf6 attr
904  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
905  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
906  end do
907 
908  ! exit define mode and write data
909  call nf_verify(nf90_enddef(ncid), nc_fname)
910  do k = 1, disv%nlay
911  call nf_verify(nf90_put_var(ncid, var_id(k), p_mem(:, k)), nc_fname)
912  end do
913 
914  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 tdismodule, only: kper
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
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  call nf_verify(nf90_put_var(ncid, &
601  var_ids%export(1), p_mem, &
602  start=(/1, kper/), &
603  count=(/disv%ncpl, 1/)), nc_fname)
604  end if
605 
606  else
607 
608  int2d(1:disv%ncpl, 1:disv%nlay) => p_mem(1:disv%nodesuser)
609 
610  if (iper == 0) then
611  allocate (var_id(disv%nlay))
612 
613  ! reenter define mode and create variable
614  call nf_verify(nf90_redef(ncid), nc_fname)
615  do k = 1, disv%nlay
616  ! set names
617  varname = export_varname(pkgname, idt%tagname, mempath, layer=k)
618  longname = export_longname(idt%longname, pkgname, idt%tagname, &
619  mempath, layer=k)
620 
621  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
622  (/dim_ids%nmesh_face/), var_id(k)), &
623  nc_fname)
624 
625  ! apply chunking parameters
626  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
627  ! deflate and shuffle
628  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
629 
630  ! put attr
631  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
632  (/nf90_fill_int/)), nc_fname)
633  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
634  longname), nc_fname)
635 
636  ! add grid mapping and mf6 attr
637  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
638  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
639  end do
640 
641  ! exit define mode and write data
642  call nf_verify(nf90_enddef(ncid), nc_fname)
643  do k = 1, disv%nlay
644  call nf_verify(nf90_put_var(ncid, var_id(k), int2d(:, k)), nc_fname)
645  end do
646 
647  ! cleanup
648  deallocate (var_id)
649  else
650  ! timeseries, add period data
651  do k = 1, disv%nlay
652  int1d(1:disv%ncpl) => int2d(:, k)
653  call nf_verify(nf90_put_var(ncid, &
654  var_ids%export(k), int1d, &
655  start=(/1, kper/), &
656  count=(/disv%ncpl, 1/)), nc_fname)
657  end do
658  end if
659  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 664 of file DisvNCMesh.f90.

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