MODFLOW 6  version 6.7.0.dev1
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.

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

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

977  integer(I4B), intent(in) :: ncid
978  type(MeshNCDimIdType), intent(inout) :: dim_ids
979  type(MeshNCVarIdType), intent(inout) :: var_ids
980  type(DisType), pointer, intent(in) :: dis
981  real(DP), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
982  character(len=*), intent(in) :: nc_varname
983  character(len=*), intent(in) :: pkgname
984  character(len=*), intent(in) :: tagname
985  character(len=*), intent(in) :: gridmap_name
986  character(len=*), intent(in) :: shapestr
987  character(len=*), intent(in) :: longname
988  character(len=*), intent(in) :: nc_tag
989  integer(I4B), intent(in) :: deflate
990  integer(I4B), intent(in) :: shuffle
991  integer(I4B), intent(in) :: chunk_face
992  character(len=*), intent(in) :: nc_fname
993  integer(I4B) :: var_id
994  character(len=LINELENGTH) :: longname_l, varname_l
995  real(DP), dimension(:), pointer, contiguous :: dbl1d
996  integer(I4B), dimension(1) :: layer_shape
997 
998  ! set names
999  varname_l = export_varname(nc_varname)
1000  longname_l = export_longname(longname, pkgname, tagname, 0)
1001 
1002  ! reenter define mode and create variable
1003  call nf_verify(nf90_redef(ncid), nc_fname)
1004  call nf_verify(nf90_def_var(ncid, varname_l, nf90_double, &
1005  (/dim_ids%nmesh_face/), var_id), &
1006  nc_fname)
1007 
1008  ! apply chunking parameters
1009  call ncvar_chunk(ncid, var_id, chunk_face, nc_fname)
1010  ! deflate and shuffle
1011  call ncvar_deflate(ncid, var_id, deflate, shuffle, nc_fname)
1012 
1013  ! put attr
1014  call nf_verify(nf90_put_att(ncid, var_id, '_FillValue', &
1015  (/nf90_fill_double/)), nc_fname)
1016  call nf_verify(nf90_put_att(ncid, var_id, 'long_name', &
1017  longname_l), nc_fname)
1018 
1019  ! add grid mapping and mf6 attr
1020  call ncvar_gridmap(ncid, var_id, gridmap_name, nc_fname)
1021  call ncvar_mf6attr(ncid, var_id, 0, 0, 0, nc_tag, nc_fname)
1022 
1023  ! exit define mode and write data
1024  call nf_verify(nf90_enddef(ncid), nc_fname)
1025  layer_shape(1) = dis%nrow * dis%ncol
1026  dbl1d(1:layer_shape(1)) => p_mem
1027  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 1032 of file DisNCMesh.f90.

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

Definition at line 598 of file DisNCMesh.f90.

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

727  integer(I4B), intent(in) :: ncid
728  type(MeshNCDimIdType), intent(inout) :: dim_ids
729  type(MeshNCVarIdType), intent(inout) :: var_ids
730  type(DisType), pointer, intent(in) :: dis
731  integer(I4B), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
732  character(len=*), intent(in) :: nc_varname
733  character(len=*), intent(in) :: pkgname
734  character(len=*), intent(in) :: tagname
735  character(len=*), intent(in) :: gridmap_name
736  character(len=*), intent(in) :: shapestr
737  character(len=*), intent(in) :: longname
738  character(len=*), intent(in) :: nc_tag
739  integer(I4B), intent(in) :: deflate
740  integer(I4B), intent(in) :: shuffle
741  integer(I4B), intent(in) :: chunk_face
742  character(len=*), intent(in) :: nc_fname
743  integer(I4B) :: var_id
744  integer(I4B), dimension(:), pointer, contiguous :: int1d
745  integer(I4B), dimension(1) :: layer_shape
746  character(len=LINELENGTH) :: longname_l, varname_l
747 
748  ! set names
749  varname_l = export_varname(nc_varname)
750  longname_l = export_longname(longname, pkgname, tagname, 0)
751 
752  ! reenter define mode and create variable
753  call nf_verify(nf90_redef(ncid), nc_fname)
754  call nf_verify(nf90_def_var(ncid, varname_l, nf90_int, &
755  (/dim_ids%nmesh_face/), var_id), &
756  nc_fname)
757 
758  ! apply chunking parameters
759  call ncvar_chunk(ncid, var_id, chunk_face, nc_fname)
760  ! deflate and shuffle
761  call ncvar_deflate(ncid, var_id, deflate, shuffle, nc_fname)
762 
763  ! put attr
764  call nf_verify(nf90_put_att(ncid, var_id, '_FillValue', &
765  (/nf90_fill_int/)), nc_fname)
766  call nf_verify(nf90_put_att(ncid, var_id, 'long_name', &
767  longname_l), nc_fname)
768 
769  ! add grid mapping and mf6 attr
770  call ncvar_gridmap(ncid, var_id, gridmap_name, nc_fname)
771  call ncvar_mf6attr(ncid, var_id, 0, 0, 0, nc_tag, nc_fname)
772 
773  ! exit define mode and write data
774  call nf_verify(nf90_enddef(ncid), nc_fname)
775  layer_shape(1) = dis%nrow * dis%ncol
776  int1d(1:layer_shape(1)) => p_mem
777  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 782 of file DisNCMesh.f90.

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