MODFLOW 6  version 6.7.0.dev3
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  ! write face nodes array to netcdf file
501  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_nodes, &
502  verts, start=(/1, n/), &
503  count=(/maxvert, 1/)), &
504  this%nc_fname)
505 
506  ! set face y bounds array
507  bnds = nf90_fill_double
508  do m = 1, size(bnds)
509  if (verts(m) /= nf90_fill_int) then
510  bnds(m) = vert_yt(verts(m))
511  end if
512  ! write face y bounds array to netcdf file
513  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_ybnds, &
514  bnds, start=(/1, n/), &
515  count=(/maxvert, 1/)), &
516  this%nc_fname)
517  end do
518 
519  ! set face x bounds array
520  bnds = nf90_fill_double
521  do m = 1, size(bnds)
522  if (verts(m) /= nf90_fill_int) then
523  bnds(m) = vert_xt(verts(m))
524  end if
525  ! write face x bounds array to netcdf file
526  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_xbnds, &
527  bnds, start=(/1, n/), &
528  count=(/maxvert, 1/)), &
529  this%nc_fname)
530  end do
531  end do
532 
533  ! cleanup
534  deallocate (bnds)
535  deallocate (verts)
536  deallocate (cell_xt)
537  deallocate (cell_yt)
538  deallocate (vert_xt)
539  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 733 of file DisvNCMesh.f90.

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

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

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

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