MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
meshdismodelmodule Module Reference

This module contains the MeshDisModelModule. More...

Data Types

type  mesh2ddisexporttype
 

Functions/Subroutines

subroutine dis_export_init (this, modelname, modeltype, modelfname, nc_fname, disenum, nctype, iout)
 netcdf export dis init More...
 
subroutine dis_export_destroy (this)
 netcdf export dis destroy More...
 
subroutine df (this)
 netcdf export define More...
 
subroutine step (this)
 netcdf export step More...
 
subroutine package_step_ilayer (this, export_pkg, ilayer_varname, ilayer)
 netcdf export package dynamic input with ilayer index variable More...
 
subroutine package_step (this, export_pkg)
 netcdf export package dynamic input More...
 
subroutine export_layer_3d (this, export_pkg, idt, ilayer_read, ialayer, dbl1d, nc_varname, input_attr, iaux)
 export layer variable as full grid 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 (ncid, dim_ids, x_dim, y_dim, var_ids, dis, p_mem, nc_varname, pkgname, tagname, gridmap_name, shapestr, longname, nc_tag, deflate, shuffle, chunk_face, iper, nc_fname)
 netcdf export 1D integer More...
 
subroutine nc_export_int2d (ncid, dim_ids, var_ids, dis, p_mem, nc_varname, pkgname, tagname, gridmap_name, shapestr, longname, nc_tag, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 2D integer More...
 
subroutine nc_export_int3d (ncid, dim_ids, var_ids, dis, p_mem, nc_varname, pkgname, tagname, gridmap_name, shapestr, longname, nc_tag, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 3D integer More...
 
subroutine nc_export_dbl1d (ncid, dim_ids, x_dim, y_dim, var_ids, dis, p_mem, nc_varname, pkgname, tagname, gridmap_name, shapestr, longname, nc_tag, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 1D double More...
 
subroutine nc_export_dbl2d (ncid, dim_ids, var_ids, dis, p_mem, nc_varname, pkgname, tagname, gridmap_name, shapestr, longname, nc_tag, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 2D double More...
 
subroutine nc_export_dbl3d (ncid, dim_ids, var_ids, dis, p_mem, nc_varname, pkgname, tagname, gridmap_name, shapestr, longname, nc_tag, deflate, shuffle, chunk_face, iper, iaux, nc_fname)
 netcdf export 3D double More...
 

Detailed Description

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

Function/Subroutine Documentation

◆ add_mesh_data()

subroutine meshdismodelmodule::add_mesh_data ( class(mesh2ddisexporttype), intent(inout)  this)

Definition at line 459 of file DisNCMesh.f90.

460  class(Mesh2dDisExportType), intent(inout) :: this
461  integer(I4B) :: cnt, maxvert, m
462  integer(I4B), dimension(:), allocatable :: verts
463  real(DP), dimension(:), allocatable :: bnds
464  integer(I4B) :: i, j
465  real(DP) :: x, y
466  real(DP), dimension(:), allocatable :: node_x, node_y
467  real(DP), dimension(:), allocatable :: cell_x, cell_y
468 
469  ! initialize max vertices required to define cell
470  maxvert = 4
471 
472  ! set mesh container variable value to 1
473  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh, 1), &
474  this%nc_fname)
475 
476  ! allocate temporary arrays
477  allocate (verts(maxvert))
478  allocate (bnds(maxvert))
479  allocate (node_x(((this%dis%ncol + 1) * (this%dis%nrow + 1))))
480  allocate (node_y(((this%dis%ncol + 1) * (this%dis%nrow + 1))))
481  allocate (cell_x((this%dis%ncol * this%dis%nrow)))
482  allocate (cell_y((this%dis%ncol * this%dis%nrow)))
483 
484  ! set node_x and node_y arrays
485  cnt = 0
486  node_x = nf90_fill_double
487  node_y = nf90_fill_double
488  y = this%dis%yorigin + sum(this%dis%delc)
489  do j = this%dis%nrow, 0, -1
490  x = this%dis%xorigin
491  do i = this%dis%ncol, 0, -1
492  cnt = cnt + 1
493  node_x(cnt) = x
494  node_y(cnt) = y
495  if (i > 0) x = x + this%dis%delr(i)
496  end do
497  if (j > 0) y = y - this%dis%delc(j)
498  end do
499 
500  ! write node_x and node_y arrays to netcdf file
501  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_x, node_x), &
502  this%nc_fname)
503  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_y, node_y), &
504  this%nc_fname)
505 
506  ! set cell_x and cell_y arrays
507  cnt = 1
508  cell_x = nf90_fill_double
509  cell_y = nf90_fill_double
510  do j = 1, this%dis%nrow
511  x = this%dis%xorigin
512  y = this%dis%celly(j) + this%dis%yorigin
513  do i = 1, this%dis%ncol
514  cell_x(cnt) = x
515  cell_y(cnt) = y
516  cnt = cnt + 1
517  x = this%dis%cellx(i) + this%dis%xorigin
518  end do
519  end do
520 
521  ! write face_x and face_y arrays to netcdf file
522  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_x, cell_x), &
523  this%nc_fname)
524  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_y, cell_y), &
525  this%nc_fname)
526 
527  ! set face nodes array
528  cnt = 0
529  do i = 1, this%dis%nrow
530  do j = 1, this%dis%ncol
531  cnt = cnt + 1
532  verts = nf90_fill_int
533  verts(1) = cnt + this%dis%ncol + i
534  verts(2) = cnt + this%dis%ncol + i + 1
535  if (i > 1) then
536  verts(3) = cnt + i
537  verts(4) = cnt + i - 1
538  else
539  verts(3) = cnt + 1
540  verts(4) = cnt
541  end if
542 
543  ! write face nodes array to netcdf file
544  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_nodes, &
545  verts, start=(/1, cnt/), &
546  count=(/maxvert, 1/)), &
547  this%nc_fname)
548 
549  ! set face y bounds array
550  bnds = nf90_fill_double
551  do m = 1, size(bnds)
552  if (verts(m) /= nf90_fill_int) then
553  bnds(m) = node_y(verts(m))
554  end if
555  ! write face y bounds array to netcdf file
556  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_ybnds, &
557  bnds, start=(/1, cnt/), &
558  count=(/maxvert, 1/)), &
559  this%nc_fname)
560  end do
561 
562  ! set face x bounds array
563  bnds = nf90_fill_double
564  do m = 1, size(bnds)
565  if (verts(m) /= nf90_fill_int) then
566  bnds(m) = node_x(verts(m))
567  end if
568  ! write face x bounds array to netcdf file
569  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_xbnds, &
570  bnds, start=(/1, cnt/), &
571  count=(/maxvert, 1/)), &
572  this%nc_fname)
573  end do
574  end do
575  end do
576 
577  ! cleanup
578  deallocate (bnds)
579  deallocate (verts)
580  deallocate (node_x)
581  deallocate (node_y)
582  deallocate (cell_x)
583  deallocate (cell_y)
Here is the call graph for this function:

◆ define_dim()

subroutine meshdismodelmodule::define_dim ( class(mesh2ddisexporttype), intent(inout)  this)
private

Definition at line 413 of file DisNCMesh.f90.

414  use constantsmodule, only: mvalidate
415  use simvariablesmodule, only: isim_mode
416  class(Mesh2dDisExportType), intent(inout) :: this
417 
418  ! time
419  if (isim_mode /= mvalidate) then
420  call nf_verify(nf90_def_dim(this%ncid, 'time', this%totnstp, &
421  this%dim_ids%time), this%nc_fname)
422  call nf_verify(nf90_def_var(this%ncid, 'time', nf90_double, &
423  this%dim_ids%time, this%var_ids%time), &
424  this%nc_fname)
425  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'calendar', &
426  'standard'), this%nc_fname)
427  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'units', &
428  this%datetime), this%nc_fname)
429  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'axis', 'T'), &
430  this%nc_fname)
431  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'standard_name', &
432  'time'), this%nc_fname)
433  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'long_name', &
434  'time'), this%nc_fname)
435  end if
436 
437  ! mesh
438  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_node', &
439  ((this%dis%ncol + 1) * (this%dis%nrow + 1)), &
440  this%dim_ids%nmesh_node), this%nc_fname)
441  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_face', &
442  (this%dis%ncol * this%dis%nrow), &
443  this%dim_ids%nmesh_face), this%nc_fname)
444  call nf_verify(nf90_def_dim(this%ncid, 'max_nmesh_face_nodes', 4, &
445  this%dim_ids%max_nmesh_face_nodes), &
446  this%nc_fname)
447 
448  ! x, y, nlay
449  call nf_verify(nf90_def_dim(this%ncid, 'nlay', this%dis%nlay, &
450  this%dim_ids%nlay), this%nc_fname)
451  call nf_verify(nf90_def_dim(this%ncid, 'x', this%dis%ncol, &
452  this%x_dim), this%nc_fname)
453  call nf_verify(nf90_def_dim(this%ncid, 'y', this%dis%nrow, &
454  this%y_dim), 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 meshdismodelmodule::df ( class(mesh2ddisexporttype), intent(inout)  this)
private

Definition at line 87 of file DisNCMesh.f90.

88  use constantsmodule, only: mvalidate
89  use simvariablesmodule, only: isim_mode
90  class(Mesh2dDisExportType), intent(inout) :: this
91  ! put root group file scope attributes
92  call this%add_global_att()
93  ! define root group dimensions and coordinate variables
94  call this%define_dim()
95  ! define mesh variables
96  call this%create_mesh()
97  if (isim_mode /= mvalidate) then
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  ! define and set package input griddata
106  call this%add_pkg_data()
107  ! define and set gridmap variable
108  call this%define_gridmap()
109  ! synchronize file
110  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
Here is the call graph for this function:

◆ dis_export_destroy()

subroutine meshdismodelmodule::dis_export_destroy ( class(mesh2ddisexporttype), intent(inout)  this)

Definition at line 77 of file DisNCMesh.f90.

78  class(Mesh2dDisExportType), intent(inout) :: this
79  deallocate (this%var_ids%dependent)
80  ! destroy base class
81  call this%mesh_destroy()
82  call this%NCModelExportType%destroy()

◆ dis_export_init()

subroutine meshdismodelmodule::dis_export_init ( class(mesh2ddisexporttype), 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 52 of file DisNCMesh.f90.

55  class(Mesh2dDisExportType), intent(inout) :: this
56  character(len=*), intent(in) :: modelname
57  character(len=*), intent(in) :: modeltype
58  character(len=*), intent(in) :: modelfname
59  character(len=*), intent(in) :: nc_fname
60  integer(I4B), intent(in) :: disenum
61  integer(I4B), intent(in) :: nctype
62  integer(I4B), intent(in) :: iout
63 
64  ! set nlay
65  this%nlay = this%dis%nlay
66 
67  ! allocate var_id arrays
68  allocate (this%var_ids%dependent(this%nlay))
69 
70  ! initialize base class
71  call this%mesh_init(modelname, modeltype, modelfname, nc_fname, disenum, &
72  nctype, iout)

◆ export_input_array()

subroutine meshdismodelmodule::export_input_array ( class(mesh2ddisexporttype), 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 340 of file DisNCMesh.f90.

341  class(Mesh2dDisExportType), intent(inout) :: this
342  character(len=*), intent(in) :: pkgtype
343  character(len=*), intent(in) :: pkgname
344  character(len=*), intent(in) :: mempath
345  type(InputParamDefinitionType), pointer, intent(in) :: idt
346  integer(I4B), dimension(:), pointer, contiguous :: int1d
347  integer(I4B), dimension(:, :), pointer, contiguous :: int2d
348  integer(I4B), dimension(:, :, :), pointer, contiguous :: int3d
349  real(DP), dimension(:), pointer, contiguous :: dbl1d
350  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
351  real(DP), dimension(:, :, :), pointer, contiguous :: dbl3d
352  character(len=LINELENGTH) :: nc_varname, input_attr
353  integer(I4B) :: iper, iaux
354 
355  iper = 0
356  iaux = 0
357 
358  ! set package base name
359  nc_varname = trim(pkgname)//'_'//trim(idt%mf6varname)
360  ! put input attributes
361  input_attr = this%input_attribute(pkgname, idt)
362 
363  select case (idt%datatype)
364  case ('INTEGER1D')
365  call mem_setptr(int1d, idt%mf6varname, mempath)
366  call nc_export_int1d(this%ncid, this%dim_ids, this%x_dim, this%y_dim, &
367  this%var_ids, this%dis, int1d, nc_varname, pkgname, &
368  idt%tagname, this%gridmap_name, idt%shape, &
369  idt%longname, input_attr, this%deflate, this%shuffle, &
370  this%chunk_face, iper, this%nc_fname)
371  case ('INTEGER2D')
372  call mem_setptr(int2d, idt%mf6varname, mempath)
373  call nc_export_int2d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
374  int2d, nc_varname, pkgname, idt%tagname, &
375  this%gridmap_name, idt%shape, idt%longname, &
376  input_attr, this%deflate, this%shuffle, &
377  this%chunk_face, this%nc_fname)
378  case ('INTEGER3D')
379  call mem_setptr(int3d, idt%mf6varname, mempath)
380  call nc_export_int3d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
381  int3d, nc_varname, pkgname, idt%tagname, &
382  this%gridmap_name, idt%shape, idt%longname, &
383  input_attr, this%deflate, this%shuffle, &
384  this%chunk_face, this%nc_fname)
385  case ('DOUBLE1D')
386  call mem_setptr(dbl1d, idt%mf6varname, mempath)
387  call nc_export_dbl1d(this%ncid, this%dim_ids, this%x_dim, this%y_dim, &
388  this%var_ids, this%dis, dbl1d, nc_varname, pkgname, &
389  idt%tagname, this%gridmap_name, idt%shape, &
390  idt%longname, input_attr, this%deflate, this%shuffle, &
391  this%chunk_face, this%nc_fname)
392  case ('DOUBLE2D')
393  call mem_setptr(dbl2d, idt%mf6varname, mempath)
394  call nc_export_dbl2d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
395  dbl2d, nc_varname, pkgname, idt%tagname, &
396  this%gridmap_name, idt%shape, idt%longname, &
397  input_attr, this%deflate, this%shuffle, &
398  this%chunk_face, this%nc_fname)
399  case ('DOUBLE3D')
400  call mem_setptr(dbl3d, idt%mf6varname, mempath)
401  call nc_export_dbl3d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
402  dbl3d, nc_varname, pkgname, idt%tagname, &
403  this%gridmap_name, idt%shape, idt%longname, &
404  input_attr, this%deflate, this%shuffle, &
405  this%chunk_face, iper, iaux, this%nc_fname)
406  case default
407  ! no-op, no other datatypes exported
408  end select
Here is the call graph for this function:

◆ export_layer_3d()

subroutine meshdismodelmodule::export_layer_3d ( class(mesh2ddisexporttype), intent(inout)  this,
class(exportpackagetype), intent(in), pointer  export_pkg,
type(inputparamdefinitiontype), intent(in), pointer  idt,
logical(lgp), intent(in)  ilayer_read,
integer(i4b), dimension(:), intent(in), pointer, contiguous  ialayer,
real(dp), dimension(:), intent(in), pointer, contiguous  dbl1d,
character(len=*), intent(in)  nc_varname,
character(len=*), intent(in)  input_attr,
integer(i4b), intent(in), optional  iaux 
)

Definition at line 282 of file DisNCMesh.f90.

284  use constantsmodule, only: dnodata, dzero
286  class(Mesh2dDisExportType), intent(inout) :: this
287  class(ExportPackageType), pointer, intent(in) :: export_pkg
288  type(InputParamDefinitionType), pointer, intent(in) :: idt
289  logical(LGP), intent(in) :: ilayer_read
290  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: ialayer
291  real(DP), dimension(:), pointer, contiguous, intent(in) :: dbl1d
292  character(len=*), intent(in) :: nc_varname
293  character(len=*), intent(in) :: input_attr
294  integer(I4B), optional, intent(in) :: iaux
295  real(DP), dimension(:, :, :), pointer, contiguous :: dbl3d
296  integer(I4B) :: n, i, j, k, nvals, idxaux
297  real(DP), dimension(:, :), contiguous, pointer :: dbl2d_ptr
298 
299  ! initialize
300  idxaux = 0
301  if (present(iaux)) then
302  idxaux = iaux
303  end if
304 
305  allocate (dbl3d(export_pkg%mshape(3), export_pkg%mshape(2), &
306  export_pkg%mshape(1)))
307 
308  if (ilayer_read) then
309  do k = 1, size(dbl3d, dim=3)
310  n = 0
311  do i = 1, size(dbl3d, dim=2)
312  do j = 1, size(dbl3d, dim=1)
313  n = n + 1
314  if (ialayer(n) == k) then
315  dbl3d(j, i, k) = dbl1d(n)
316  else
317  dbl3d(j, i, k) = dnodata
318  end if
319  end do
320  end do
321  end do
322  else
323  dbl3d = dnodata
324  nvals = export_pkg%mshape(3) * export_pkg%mshape(2)
325  dbl2d_ptr(1:export_pkg%mshape(3), 1:export_pkg%mshape(2)) => dbl1d(1:nvals)
326  dbl3d(:, :, 1) = dbl2d_ptr(:, :)
327  end if
328 
329  call nc_export_dbl3d(this%ncid, this%dim_ids, this%var_ids, this%dis, dbl3d, &
330  nc_varname, export_pkg%mf6_input%subcomponent_name, &
331  idt%tagname, this%gridmap_name, idt%shape, &
332  idt%longname, input_attr, this%deflate, this%shuffle, &
333  this%chunk_face, export_pkg%iper, idxaux, this%nc_fname)
334 
335  deallocate (dbl3d)
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
This module contains the NCModelExportModule.
Definition: NCModel.f90:8
Here is the call graph for this function:

◆ nc_export_dbl1d()

subroutine meshdismodelmodule::nc_export_dbl1d ( integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
integer(i4b), intent(in)  x_dim,
integer(i4b), intent(in)  y_dim,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
real(dp), dimension(:), intent(in), pointer, contiguous  p_mem,
character(len=*), intent(in)  nc_varname,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  tagname,
character(len=*), intent(in)  gridmap_name,
character(len=*), intent(in)  shapestr,
character(len=*), intent(in)  longname,
character(len=*), intent(in)  nc_tag,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)
private

Definition at line 840 of file DisNCMesh.f90.

844  integer(I4B), intent(in) :: ncid
845  type(MeshNCDimIdType), intent(inout) :: dim_ids
846  integer(I4B), intent(in) :: x_dim
847  integer(I4B), intent(in) :: y_dim
848  type(MeshNCVarIdType), intent(inout) :: var_ids
849  type(DisType), pointer, intent(in) :: dis
850  real(DP), dimension(:), pointer, contiguous, intent(in) :: p_mem
851  character(len=*), intent(in) :: nc_varname
852  character(len=*), intent(in) :: pkgname
853  character(len=*), intent(in) :: tagname
854  character(len=*), intent(in) :: gridmap_name
855  character(len=*), intent(in) :: shapestr
856  character(len=*), intent(in) :: longname
857  character(len=*), intent(in) :: nc_tag
858  integer(I4B), intent(in) :: deflate
859  integer(I4B), intent(in) :: shuffle
860  integer(I4B), intent(in) :: chunk_face
861  character(len=*), intent(in) :: nc_fname
862  integer(I4B), dimension(3) :: dis_shape
863  integer(I4B), dimension(1) :: layer_shape
864  real(DP), dimension(:, :, :), pointer, contiguous :: dbl3d
865  real(DP), dimension(:), pointer, contiguous :: dbl1d
866  integer(I4B) :: axis_dim, nvals, k
867  integer(I4B), dimension(:), allocatable :: var_id
868  character(len=LINELENGTH) :: longname_l, varname_l
869 
870  if (shapestr == 'NROW' .or. &
871  shapestr == 'NCOL') then ! .or. &
872  !shapestr == 'NCPL') then
873 
874  select case (shapestr)
875  case ('NROW')
876  axis_dim = y_dim
877  case ('NCOL')
878  axis_dim = x_dim
879  !case ('NCPL')
880  ! axis_dim = dim_ids%nmesh_face
881  end select
882 
883  ! set names
884  varname_l = export_varname(nc_varname)
885  longname_l = export_longname(longname, pkgname, tagname, 0)
886 
887  allocate (var_id(1))
888 
889  ! reenter define mode and create variable
890  call nf_verify(nf90_redef(ncid), nc_fname)
891  call nf_verify(nf90_def_var(ncid, varname_l, nf90_double, &
892  (/axis_dim/), var_id(1)), &
893  nc_fname)
894 
895  ! NROW/NCOL shapes use default chunking
896  call ncvar_deflate(ncid, var_id(1), deflate, shuffle, nc_fname)
897 
898  ! put attr
899  call nf_verify(nf90_put_att(ncid, var_id(1), '_FillValue', &
900  (/nf90_fill_double/)), nc_fname)
901  call nf_verify(nf90_put_att(ncid, var_id(1), 'long_name', &
902  longname_l), nc_fname)
903 
904  ! add mf6 attr
905  call ncvar_mf6attr(ncid, var_id(1), 0, 0, 0, nc_tag, nc_fname)
906 
907  ! exit define mode and write data
908  call nf_verify(nf90_enddef(ncid), nc_fname)
909  call nf_verify(nf90_put_var(ncid, var_id(1), p_mem), &
910  nc_fname)
911 
912  else
913  allocate (var_id(dis%nlay))
914 
915  ! reenter define mode and create variable
916  call nf_verify(nf90_redef(ncid), nc_fname)
917  do k = 1, dis%nlay
918  ! set names
919  varname_l = export_varname(nc_varname, layer=k)
920  longname_l = export_longname(longname, pkgname, tagname, k)
921 
922  call nf_verify(nf90_def_var(ncid, varname_l, nf90_double, &
923  (/dim_ids%nmesh_face/), var_id(k)), &
924  nc_fname)
925 
926  ! apply chunking parameters
927  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
928  ! defalte and shuffle
929  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
930 
931  ! put attr
932  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
933  (/nf90_fill_double/)), nc_fname)
934  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
935  longname_l), nc_fname)
936 
937  ! add grid mapping and mf6 attr
938  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
939  call ncvar_mf6attr(ncid, var_id(k), k, 0, 0, nc_tag, nc_fname)
940  end do
941 
942  ! reshape input
943  dis_shape(1) = dis%ncol
944  dis_shape(2) = dis%nrow
945  dis_shape(3) = dis%nlay
946  nvals = product(dis_shape)
947  dbl3d(1:dis_shape(1), 1:dis_shape(2), 1:dis_shape(3)) => p_mem(1:nvals)
948 
949  ! exit define mode and write data
950  call nf_verify(nf90_enddef(ncid), nc_fname)
951  layer_shape(1) = dis%nrow * dis%ncol
952  do k = 1, dis%nlay
953  dbl1d(1:layer_shape(1)) => dbl3d(:, :, k)
954  call nf_verify(nf90_put_var(ncid, var_id(k), dbl1d), nc_fname)
955  end do
956 
957  ! cleanup
958  deallocate (var_id)
959  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_dbl2d()

subroutine meshdismodelmodule::nc_export_dbl2d ( integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
real(dp), dimension(:, :), intent(in), pointer, contiguous  p_mem,
character(len=*), intent(in)  nc_varname,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  tagname,
character(len=*), intent(in)  gridmap_name,
character(len=*), intent(in)  shapestr,
character(len=*), intent(in)  longname,
character(len=*), intent(in)  nc_tag,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)
private

Definition at line 964 of file DisNCMesh.f90.

967  integer(I4B), intent(in) :: ncid
968  type(MeshNCDimIdType), intent(inout) :: dim_ids
969  type(MeshNCVarIdType), intent(inout) :: var_ids
970  type(DisType), pointer, intent(in) :: dis
971  real(DP), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
972  character(len=*), intent(in) :: nc_varname
973  character(len=*), intent(in) :: pkgname
974  character(len=*), intent(in) :: tagname
975  character(len=*), intent(in) :: gridmap_name
976  character(len=*), intent(in) :: shapestr
977  character(len=*), intent(in) :: longname
978  character(len=*), intent(in) :: nc_tag
979  integer(I4B), intent(in) :: deflate
980  integer(I4B), intent(in) :: shuffle
981  integer(I4B), intent(in) :: chunk_face
982  character(len=*), intent(in) :: nc_fname
983  integer(I4B) :: var_id
984  character(len=LINELENGTH) :: longname_l, varname_l
985  real(DP), dimension(:), pointer, contiguous :: dbl1d
986  integer(I4B), dimension(1) :: layer_shape
987 
988  ! set names
989  varname_l = export_varname(nc_varname)
990  longname_l = export_longname(longname, pkgname, tagname, 0)
991 
992  ! reenter define mode and create variable
993  call nf_verify(nf90_redef(ncid), nc_fname)
994  call nf_verify(nf90_def_var(ncid, varname_l, nf90_double, &
995  (/dim_ids%nmesh_face/), var_id), &
996  nc_fname)
997 
998  ! apply chunking parameters
999  call ncvar_chunk(ncid, var_id, chunk_face, nc_fname)
1000  ! deflate and shuffle
1001  call ncvar_deflate(ncid, var_id, deflate, shuffle, nc_fname)
1002 
1003  ! put attr
1004  call nf_verify(nf90_put_att(ncid, var_id, '_FillValue', &
1005  (/nf90_fill_double/)), nc_fname)
1006  call nf_verify(nf90_put_att(ncid, var_id, 'long_name', &
1007  longname_l), nc_fname)
1008 
1009  ! add grid mapping and mf6 attr
1010  call ncvar_gridmap(ncid, var_id, gridmap_name, nc_fname)
1011  call ncvar_mf6attr(ncid, var_id, 0, 0, 0, nc_tag, nc_fname)
1012 
1013  ! exit define mode and write data
1014  call nf_verify(nf90_enddef(ncid), nc_fname)
1015  layer_shape(1) = dis%nrow * dis%ncol
1016  dbl1d(1:layer_shape(1)) => p_mem
1017  call nf_verify(nf90_put_var(ncid, var_id, dbl1d), nc_fname)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_dbl3d()

subroutine meshdismodelmodule::nc_export_dbl3d ( integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
real(dp), dimension(:, :, :), intent(in), pointer, contiguous  p_mem,
character(len=*), intent(in)  nc_varname,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  tagname,
character(len=*), intent(in)  gridmap_name,
character(len=*), intent(in)  shapestr,
character(len=*), intent(in)  longname,
character(len=*), intent(in)  nc_tag,
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 1022 of file DisNCMesh.f90.

1026  use constantsmodule, only: dnodata
1027  integer(I4B), intent(in) :: ncid
1028  type(MeshNCDimIdType), intent(inout) :: dim_ids
1029  type(MeshNCVarIdType), intent(inout) :: var_ids
1030  type(DisType), pointer, intent(in) :: dis
1031  real(DP), dimension(:, :, :), pointer, contiguous, intent(in) :: p_mem
1032  character(len=*), intent(in) :: nc_varname
1033  character(len=*), intent(in) :: pkgname
1034  character(len=*), intent(in) :: tagname
1035  character(len=*), intent(in) :: gridmap_name
1036  character(len=*), intent(in) :: shapestr
1037  character(len=*), intent(in) :: longname
1038  character(len=*), intent(in) :: nc_tag
1039  integer(I4B), intent(in) :: deflate
1040  integer(I4B), intent(in) :: shuffle
1041  integer(I4B), intent(in) :: chunk_face
1042  integer(I4B), intent(in) :: iper
1043  integer(I4B), intent(in) :: iaux
1044  character(len=*), intent(in) :: nc_fname
1045  integer(I4B), dimension(:), allocatable :: var_id
1046  real(DP), dimension(:), pointer, contiguous :: dbl1d
1047  character(len=LINELENGTH) :: longname_l, varname_l
1048  integer(I4B), dimension(1) :: layer_shape
1049  integer(I4B) :: k
1050  real(DP) :: fill_value
1051 
1052  if (iper > 0) then
1053  fill_value = dnodata
1054  else
1055  fill_value = nf90_fill_double
1056  end if
1057 
1058  allocate (var_id(dis%nlay))
1059 
1060  ! reenter define mode and create variable
1061  call nf_verify(nf90_redef(ncid), nc_fname)
1062  do k = 1, dis%nlay
1063  ! set names
1064  varname_l = export_varname(nc_varname, layer=k, iper=iper, iaux=iaux)
1065  longname_l = export_longname(longname, pkgname, tagname, layer=k, iper=iper)
1066 
1067  call nf_verify(nf90_def_var(ncid, varname_l, nf90_double, &
1068  (/dim_ids%nmesh_face/), var_id(k)), &
1069  nc_fname)
1070 
1071  ! apply chunking parameters
1072  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
1073  ! deflate and shuffle
1074  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
1075 
1076  ! put attr
1077  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
1078  (/fill_value/)), nc_fname)
1079  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
1080  longname_l), nc_fname)
1081 
1082  ! add grid mapping and mf6 attr
1083  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
1084  call ncvar_mf6attr(ncid, var_id(k), k, iper, iaux, nc_tag, nc_fname)
1085  !end if
1086  end do
1087 
1088  ! exit define mode and write data
1089  call nf_verify(nf90_enddef(ncid), nc_fname)
1090  layer_shape(1) = dis%nrow * dis%ncol
1091  do k = 1, dis%nlay
1092  dbl1d(1:layer_shape(1)) => p_mem(:, :, k)
1093  call nf_verify(nf90_put_var(ncid, var_id(k), dbl1d), nc_fname)
1094  end do
1095 
1096  ! cleanup
1097  deallocate (var_id)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_int1d()

subroutine meshdismodelmodule::nc_export_int1d ( integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
integer(i4b), intent(in)  x_dim,
integer(i4b), intent(in)  y_dim,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
integer(i4b), dimension(:), intent(in), pointer, contiguous  p_mem,
character(len=*), intent(in)  nc_varname,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  tagname,
character(len=*), intent(in)  gridmap_name,
character(len=*), intent(in)  shapestr,
character(len=*), intent(in)  longname,
character(len=*), intent(in)  nc_tag,
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 
)
private

Definition at line 588 of file DisNCMesh.f90.

592  integer(I4B), intent(in) :: ncid
593  type(MeshNCDimIdType), intent(inout) :: dim_ids
594  integer(I4B), intent(in) :: x_dim
595  integer(I4B), intent(in) :: y_dim
596  type(MeshNCVarIdType), intent(inout) :: var_ids
597  type(DisType), pointer, intent(in) :: dis
598  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: p_mem
599  character(len=*), intent(in) :: nc_varname
600  character(len=*), intent(in) :: pkgname
601  character(len=*), intent(in) :: tagname
602  character(len=*), intent(in) :: gridmap_name
603  character(len=*), intent(in) :: shapestr
604  character(len=*), intent(in) :: longname
605  character(len=*), intent(in) :: nc_tag
606  integer(I4B), intent(in) :: deflate
607  integer(I4B), intent(in) :: shuffle
608  integer(I4B), intent(in) :: chunk_face
609  integer(I4B), intent(in) :: iper
610  character(len=*), intent(in) :: nc_fname
611  integer(I4B), dimension(3) :: dis_shape
612  integer(I4B), dimension(1) :: layer_shape
613  integer(I4B), dimension(:, :, :), pointer, contiguous :: int3d
614  integer(I4B), dimension(:), pointer, contiguous :: int1d
615  integer(I4B) :: axis_dim, nvals, k
616  integer(I4B), dimension(:), allocatable :: var_id
617  character(len=LINELENGTH) :: longname_l, varname_l
618 
619  if (shapestr == 'NROW' .or. &
620  shapestr == 'NCOL' .or. &
621  shapestr == 'NCPL') then
622 
623  select case (shapestr)
624  case ('NROW')
625  axis_dim = y_dim
626  case ('NCOL')
627  axis_dim = x_dim
628  case ('NCPL')
629  axis_dim = dim_ids%nmesh_face
630  end select
631 
632  ! set names
633  varname_l = export_varname(nc_varname, layer=0, iper=iper)
634  longname_l = export_longname(longname, pkgname, tagname, layer=0, iper=iper)
635 
636  allocate (var_id(1))
637 
638  ! reenter define mode and create variable
639  call nf_verify(nf90_redef(ncid), nc_fname)
640  call nf_verify(nf90_def_var(ncid, varname_l, nf90_int, &
641  (/axis_dim/), var_id(1)), &
642  nc_fname)
643 
644  ! NROW/NCOL shapes use default chunking
645  call ncvar_deflate(ncid, var_id(1), deflate, shuffle, nc_fname)
646 
647  ! put attr
648  call nf_verify(nf90_put_att(ncid, var_id(1), '_FillValue', &
649  (/nf90_fill_int/)), nc_fname)
650  call nf_verify(nf90_put_att(ncid, var_id(1), 'long_name', &
651  longname_l), nc_fname)
652 
653  ! add mf6 attr
654  call ncvar_mf6attr(ncid, var_id(1), 0, iper, 0, nc_tag, nc_fname)
655 
656  ! exit define mode and write data
657  call nf_verify(nf90_enddef(ncid), nc_fname)
658  call nf_verify(nf90_put_var(ncid, var_id(1), p_mem), &
659  nc_fname)
660 
661  else
662  allocate (var_id(dis%nlay))
663 
664  ! reenter define mode and create variable
665  call nf_verify(nf90_redef(ncid), nc_fname)
666  do k = 1, dis%nlay
667  ! set names
668  varname_l = export_varname(nc_varname, layer=k, iper=iper)
669  longname_l = export_longname(longname, pkgname, tagname, layer=k, &
670  iper=iper)
671 
672  call nf_verify(nf90_def_var(ncid, varname_l, nf90_int, &
673  (/dim_ids%nmesh_face/), var_id(k)), &
674  nc_fname)
675 
676  ! apply chunking parameters
677  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
678  ! defalte and shuffle
679  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
680 
681  ! put attr
682  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
683  (/nf90_fill_int/)), nc_fname)
684  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
685  longname_l), nc_fname)
686 
687  ! add grid mapping and mf6 attr
688  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
689  call ncvar_mf6attr(ncid, var_id(k), k, iper, 0, nc_tag, nc_fname)
690  end do
691 
692  ! reshape input
693  dis_shape(1) = dis%ncol
694  dis_shape(2) = dis%nrow
695  dis_shape(3) = dis%nlay
696  nvals = product(dis_shape)
697  int3d(1:dis_shape(1), 1:dis_shape(2), 1:dis_shape(3)) => p_mem(1:nvals)
698 
699  ! exit define mode and write data
700  call nf_verify(nf90_enddef(ncid), nc_fname)
701  layer_shape(1) = dis%nrow * dis%ncol
702  do k = 1, dis%nlay
703  int1d(1:layer_shape(1)) => int3d(:, :, k)
704  call nf_verify(nf90_put_var(ncid, var_id(k), int1d), nc_fname)
705  end do
706 
707  ! cleanup
708  deallocate (var_id)
709  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_int2d()

subroutine meshdismodelmodule::nc_export_int2d ( integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
integer(i4b), dimension(:, :), intent(in), pointer, contiguous  p_mem,
character(len=*), intent(in)  nc_varname,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  tagname,
character(len=*), intent(in)  gridmap_name,
character(len=*), intent(in)  shapestr,
character(len=*), intent(in)  longname,
character(len=*), intent(in)  nc_tag,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)
private

Definition at line 714 of file DisNCMesh.f90.

717  integer(I4B), intent(in) :: ncid
718  type(MeshNCDimIdType), intent(inout) :: dim_ids
719  type(MeshNCVarIdType), intent(inout) :: var_ids
720  type(DisType), pointer, intent(in) :: dis
721  integer(I4B), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
722  character(len=*), intent(in) :: nc_varname
723  character(len=*), intent(in) :: pkgname
724  character(len=*), intent(in) :: tagname
725  character(len=*), intent(in) :: gridmap_name
726  character(len=*), intent(in) :: shapestr
727  character(len=*), intent(in) :: longname
728  character(len=*), intent(in) :: nc_tag
729  integer(I4B), intent(in) :: deflate
730  integer(I4B), intent(in) :: shuffle
731  integer(I4B), intent(in) :: chunk_face
732  character(len=*), intent(in) :: nc_fname
733  integer(I4B) :: var_id
734  integer(I4B), dimension(:), pointer, contiguous :: int1d
735  integer(I4B), dimension(1) :: layer_shape
736  character(len=LINELENGTH) :: longname_l, varname_l
737 
738  ! set names
739  varname_l = export_varname(nc_varname)
740  longname_l = export_longname(longname, pkgname, tagname, 0)
741 
742  ! reenter define mode and create variable
743  call nf_verify(nf90_redef(ncid), nc_fname)
744  call nf_verify(nf90_def_var(ncid, varname_l, nf90_int, &
745  (/dim_ids%nmesh_face/), var_id), &
746  nc_fname)
747 
748  ! apply chunking parameters
749  call ncvar_chunk(ncid, var_id, chunk_face, nc_fname)
750  ! deflate and shuffle
751  call ncvar_deflate(ncid, var_id, deflate, shuffle, nc_fname)
752 
753  ! put attr
754  call nf_verify(nf90_put_att(ncid, var_id, '_FillValue', &
755  (/nf90_fill_int/)), nc_fname)
756  call nf_verify(nf90_put_att(ncid, var_id, 'long_name', &
757  longname_l), nc_fname)
758 
759  ! add grid mapping and mf6 attr
760  call ncvar_gridmap(ncid, var_id, gridmap_name, nc_fname)
761  call ncvar_mf6attr(ncid, var_id, 0, 0, 0, nc_tag, nc_fname)
762 
763  ! exit define mode and write data
764  call nf_verify(nf90_enddef(ncid), nc_fname)
765  layer_shape(1) = dis%nrow * dis%ncol
766  int1d(1:layer_shape(1)) => p_mem
767  call nf_verify(nf90_put_var(ncid, var_id, int1d), nc_fname)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_int3d()

subroutine meshdismodelmodule::nc_export_int3d ( integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
integer(i4b), dimension(:, :, :), intent(in), pointer, contiguous  p_mem,
character(len=*), intent(in)  nc_varname,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  tagname,
character(len=*), intent(in)  gridmap_name,
character(len=*), intent(in)  shapestr,
character(len=*), intent(in)  longname,
character(len=*), intent(in)  nc_tag,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)
private

Definition at line 772 of file DisNCMesh.f90.

775  integer(I4B), intent(in) :: ncid
776  type(MeshNCDimIdType), intent(inout) :: dim_ids
777  type(MeshNCVarIdType), intent(inout) :: var_ids
778  type(DisType), pointer, intent(in) :: dis
779  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(in) :: p_mem
780  character(len=*), intent(in) :: nc_varname
781  character(len=*), intent(in) :: pkgname
782  character(len=*), intent(in) :: tagname
783  character(len=*), intent(in) :: gridmap_name
784  character(len=*), intent(in) :: shapestr
785  character(len=*), intent(in) :: longname
786  character(len=*), intent(in) :: nc_tag
787  integer(I4B), intent(in) :: deflate
788  integer(I4B), intent(in) :: shuffle
789  integer(I4B), intent(in) :: chunk_face
790  character(len=*), intent(in) :: nc_fname
791  integer(I4B), dimension(:), allocatable :: var_id
792  integer(I4B), dimension(:), pointer, contiguous :: int1d
793  character(len=LINELENGTH) :: longname_l, varname_l
794  integer(I4B), dimension(1) :: layer_shape
795  integer(I4B) :: k
796 
797  allocate (var_id(dis%nlay))
798 
799  ! reenter define mode and create variable
800  call nf_verify(nf90_redef(ncid), nc_fname)
801  do k = 1, dis%nlay
802  ! set names
803  varname_l = export_varname(nc_varname, layer=k)
804  longname_l = export_longname(longname, pkgname, tagname, k)
805 
806  call nf_verify(nf90_def_var(ncid, varname_l, nf90_int, &
807  (/dim_ids%nmesh_face/), var_id(k)), &
808  nc_fname)
809 
810  ! apply chunking parameters
811  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
812  ! deflate and shuffle
813  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
814 
815  ! put attr
816  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
817  (/nf90_fill_int/)), nc_fname)
818  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
819  longname_l), nc_fname)
820 
821  ! add grid mapping and mf6 attr
822  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
823  call ncvar_mf6attr(ncid, var_id(k), k, 0, 0, nc_tag, nc_fname)
824  end do
825 
826  ! exit define mode and write data
827  call nf_verify(nf90_enddef(ncid), nc_fname)
828  layer_shape(1) = dis%nrow * dis%ncol
829  do k = 1, dis%nlay
830  int1d(1:layer_shape(1)) => p_mem(:, :, k)
831  call nf_verify(nf90_put_var(ncid, var_id(k), int1d), nc_fname)
832  end do
833 
834  ! cleanup
835  deallocate (var_id)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ package_step()

subroutine meshdismodelmodule::package_step ( class(mesh2ddisexporttype), intent(inout)  this,
class(exportpackagetype), intent(in), pointer  export_pkg 
)

Definition at line 267 of file DisNCMesh.f90.

269  class(Mesh2dDisExportType), intent(inout) :: this
270  class(ExportPackageType), pointer, intent(in) :: export_pkg
271  errmsg = 'NetCDF period export not supported for model='// &
272  trim(this%modelname)//', package='// &
273  trim(export_pkg%mf6_input%subcomponent_name)
274  call store_error(errmsg, .true.)
275 
276  ! synchronize file
277  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
Here is the call graph for this function:

◆ package_step_ilayer()

subroutine meshdismodelmodule::package_step_ilayer ( class(mesh2ddisexporttype), intent(inout)  this,
class(exportpackagetype), intent(in), pointer  export_pkg,
character(len=*), intent(in)  ilayer_varname,
integer(i4b), intent(in)  ilayer 
)

Definition at line 181 of file DisNCMesh.f90.

182  use constantsmodule, only: dnodata, dzero
183  use tdismodule, only: kper
186  class(Mesh2dDisExportType), intent(inout) :: this
187  class(ExportPackageType), pointer, intent(in) :: export_pkg
188  character(len=*), intent(in) :: ilayer_varname
189  integer(I4B), intent(in) :: ilayer
190  type(InputParamDefinitionType), pointer :: idt
191  integer(I4B), dimension(:), pointer, contiguous :: int1d
192  real(DP), dimension(:), pointer, contiguous :: dbl1d
193  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
194  integer(I4B), dimension(:), pointer, contiguous :: ialayer
195  real(DP), dimension(:), contiguous, pointer :: dbl1d_ptr
196  character(len=LINELENGTH) :: nc_varname, input_attr
197  integer(I4B) :: n, iparam, nvals
198  logical(LGP) :: ilayer_read
199 
200  ! initialize
201  nullify (ialayer)
202  ilayer_read = .false.
203 
204  ! set pointer to ilayer variable
205  call mem_setptr(ialayer, export_pkg%param_names(ilayer), &
206  export_pkg%mf6_input%mempath)
207 
208  ! check if layer index variable was read
209  if (export_pkg%param_reads(ilayer)%invar == 1) then
210  ilayer_read = .true.
211  end if
212 
213  ! export defined period input
214  do iparam = 1, export_pkg%nparam
215  ! check if variable was read this period
216  if (export_pkg%param_reads(iparam)%invar < 1) cycle
217 
218  ! set input definition
219  idt => &
220  get_param_definition_type(export_pkg%mf6_input%param_dfns, &
221  export_pkg%mf6_input%component_type, &
222  export_pkg%mf6_input%subcomponent_type, &
223  'PERIOD', export_pkg%param_names(iparam), '')
224  ! set variable name and input string
225  nc_varname = trim(export_pkg%mf6_input%subcomponent_name)//'_'// &
226  trim(idt%mf6varname)
227  input_attr = this%input_attribute(export_pkg%mf6_input%subcomponent_name, &
228  idt)
229  ! export arrays
230  select case (idt%datatype)
231  case ('INTEGER1D')
232  call mem_setptr(int1d, idt%mf6varname, export_pkg%mf6_input%mempath)
233  call nc_export_int1d(this%ncid, this%dim_ids, this%x_dim, this%y_dim, &
234  this%var_ids, this%dis, int1d, nc_varname, &
235  export_pkg%mf6_input%subcomponent_name, &
236  idt%tagname, this%gridmap_name, idt%shape, &
237  idt%longname, input_attr, this%deflate, &
238  this%shuffle, this%chunk_face, kper, this%nc_fname)
239  case ('DOUBLE1D')
240  call mem_setptr(dbl1d, idt%mf6varname, export_pkg%mf6_input%mempath)
241  call this%export_layer_3d(export_pkg, idt, ilayer_read, ialayer, &
242  dbl1d, nc_varname, input_attr)
243  case ('DOUBLE2D')
244  call mem_setptr(dbl2d, idt%mf6varname, export_pkg%mf6_input%mempath)
245  nvals = this%dis%ncol * this%dis%nrow
246 
247  do n = 1, size(dbl2d, dim=1) !naux
248  dbl1d_ptr(1:nvals) => dbl2d(n, :)
249  if (all(dbl1d_ptr == dzero)) then
250  else
251  call this%export_layer_3d(export_pkg, idt, ilayer_read, ialayer, &
252  dbl1d_ptr, nc_varname, input_attr, n)
253  end if
254  end do
255  case default
256  errmsg = 'EXPORT ilayaer unsupported datatype='//trim(idt%datatype)
257  call store_error(errmsg, .true.)
258  end select
259  end do
260 
261  ! synchronize file
262  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.
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
Here is the call graph for this function:

◆ step()

subroutine meshdismodelmodule::step ( class(mesh2ddisexporttype), intent(inout)  this)

Definition at line 115 of file DisNCMesh.f90.

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