MODFLOW 6  version 6.6.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, 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 472 of file DisNCMesh.f90.

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

427  use constantsmodule, only: mvalidate
428  use simvariablesmodule, only: isim_mode
429  class(Mesh2dDisExportType), intent(inout) :: this
430  !
431  ! -- time
432  if (isim_mode /= mvalidate) then
433  call nf_verify(nf90_def_dim(this%ncid, 'time', this%totnstp, &
434  this%dim_ids%time), this%nc_fname)
435  call nf_verify(nf90_def_var(this%ncid, 'time', nf90_double, &
436  this%dim_ids%time, this%var_ids%time), &
437  this%nc_fname)
438  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'calendar', &
439  'standard'), this%nc_fname)
440  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'units', &
441  this%datetime), this%nc_fname)
442  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'axis', 'T'), &
443  this%nc_fname)
444  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'standard_name', &
445  'time'), this%nc_fname)
446  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'long_name', &
447  'time'), this%nc_fname)
448  end if
449  !
450  ! -- mesh
451  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_node', &
452  ((this%dis%ncol + 1) * (this%dis%nrow + 1)), &
453  this%dim_ids%nmesh_node), this%nc_fname)
454  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_face', &
455  (this%dis%ncol * this%dis%nrow), &
456  this%dim_ids%nmesh_face), this%nc_fname)
457  call nf_verify(nf90_def_dim(this%ncid, 'max_nmesh_face_nodes', 4, &
458  this%dim_ids%max_nmesh_face_nodes), &
459  this%nc_fname)
460  !
461  ! -- x, y, nlay
462  call nf_verify(nf90_def_dim(this%ncid, 'nlay', this%dis%nlay, &
463  this%dim_ids%nlay), this%nc_fname)
464  call nf_verify(nf90_def_dim(this%ncid, 'x', this%dis%ncol, &
465  this%x_dim), this%nc_fname)
466  call nf_verify(nf90_def_dim(this%ncid, 'y', this%dis%nrow, &
467  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 75 of file DisNCMesh.f90.

76  class(Mesh2dDisExportType), intent(inout) :: this
77  !
78  deallocate (this%var_ids%dependent)
79  !
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,
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  integer(I4B), intent(in) :: disenum
60  integer(I4B), intent(in) :: nctype
61  integer(I4B), intent(in) :: iout
62  !
63  ! -- set nlay
64  this%nlay = this%dis%nlay
65  !
66  ! allocate var_id arrays
67  allocate (this%var_ids%dependent(this%nlay))
68  !
69  ! -- initialize base class
70  call this%mesh_init(modelname, modeltype, modelfname, disenum, 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 353 of file DisNCMesh.f90.

354  class(Mesh2dDisExportType), intent(inout) :: this
355  character(len=*), intent(in) :: pkgtype
356  character(len=*), intent(in) :: pkgname
357  character(len=*), intent(in) :: mempath
358  type(InputParamDefinitionType), pointer, intent(in) :: idt
359  integer(I4B), dimension(:), pointer, contiguous :: int1d
360  integer(I4B), dimension(:, :), pointer, contiguous :: int2d
361  integer(I4B), dimension(:, :, :), pointer, contiguous :: int3d
362  real(DP), dimension(:), pointer, contiguous :: dbl1d
363  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
364  real(DP), dimension(:, :, :), pointer, contiguous :: dbl3d
365  character(len=LINELENGTH) :: nc_varname, input_attr
366  integer(I4B) :: iper, iaux
367  !
368  iper = 0
369  iaux = 0
370  !
371  ! -- set package base name
372  nc_varname = trim(pkgname)//'_'//trim(idt%mf6varname)
373  ! -- put input attributes
374  input_attr = this%input_attribute(pkgname, idt)
375  !
376  select case (idt%datatype)
377  case ('INTEGER1D')
378  call mem_setptr(int1d, idt%mf6varname, mempath)
379  call nc_export_int1d(this%ncid, this%dim_ids, this%x_dim, this%y_dim, &
380  this%var_ids, this%dis, int1d, nc_varname, pkgname, &
381  idt%tagname, this%gridmap_name, idt%shape, &
382  idt%longname, input_attr, this%deflate, this%shuffle, &
383  this%chunk_face, iper, this%nc_fname)
384  case ('INTEGER2D')
385  call mem_setptr(int2d, idt%mf6varname, mempath)
386  call nc_export_int2d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
387  int2d, nc_varname, pkgname, idt%tagname, &
388  this%gridmap_name, idt%shape, idt%longname, &
389  input_attr, this%deflate, this%shuffle, &
390  this%chunk_face, this%nc_fname)
391  case ('INTEGER3D')
392  call mem_setptr(int3d, idt%mf6varname, mempath)
393  call nc_export_int3d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
394  int3d, nc_varname, pkgname, idt%tagname, &
395  this%gridmap_name, idt%shape, idt%longname, &
396  input_attr, this%deflate, this%shuffle, &
397  this%chunk_face, this%nc_fname)
398  case ('DOUBLE1D')
399  call mem_setptr(dbl1d, idt%mf6varname, mempath)
400  call nc_export_dbl1d(this%ncid, this%dim_ids, this%x_dim, this%y_dim, &
401  this%var_ids, this%dis, dbl1d, nc_varname, pkgname, &
402  idt%tagname, this%gridmap_name, idt%shape, &
403  idt%longname, input_attr, this%deflate, this%shuffle, &
404  this%chunk_face, this%nc_fname)
405  case ('DOUBLE2D')
406  call mem_setptr(dbl2d, idt%mf6varname, mempath)
407  call nc_export_dbl2d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
408  dbl2d, nc_varname, pkgname, idt%tagname, &
409  this%gridmap_name, idt%shape, idt%longname, &
410  input_attr, this%deflate, this%shuffle, &
411  this%chunk_face, this%nc_fname)
412  case ('DOUBLE3D')
413  call mem_setptr(dbl3d, idt%mf6varname, mempath)
414  call nc_export_dbl3d(this%ncid, this%dim_ids, this%var_ids, this%dis, &
415  dbl3d, nc_varname, pkgname, idt%tagname, &
416  this%gridmap_name, idt%shape, idt%longname, &
417  input_attr, this%deflate, this%shuffle, &
418  this%chunk_face, iper, iaux, this%nc_fname)
419  case default
420  ! -- no-op, no other datatypes exported
421  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 294 of file DisNCMesh.f90.

296  use constantsmodule, only: dnodata, dzero
298  class(Mesh2dDisExportType), intent(inout) :: this
299  class(ExportPackageType), pointer, intent(in) :: export_pkg
300  type(InputParamDefinitionType), pointer, intent(in) :: idt
301  logical(LGP), intent(in) :: ilayer_read
302  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: ialayer
303  real(DP), dimension(:), pointer, contiguous, intent(in) :: dbl1d
304  character(len=*), intent(in) :: nc_varname
305  character(len=*), intent(in) :: input_attr
306  integer(I4B), optional, intent(in) :: iaux
307  ! -- local
308  real(DP), dimension(:, :, :), pointer, contiguous :: dbl3d
309  integer(I4B) :: n, i, j, k, nvals, idxaux
310  real(DP), dimension(:, :), contiguous, pointer :: dbl2d_ptr
311  !
312  ! -- initialize
313  idxaux = 0
314  if (present(iaux)) then
315  idxaux = iaux
316  end if
317 
318  allocate (dbl3d(export_pkg%mshape(3), export_pkg%mshape(2), &
319  export_pkg%mshape(1)))
320  !
321  if (ilayer_read) then
322  do k = 1, size(dbl3d, dim=3)
323  n = 0
324  do i = 1, size(dbl3d, dim=2)
325  do j = 1, size(dbl3d, dim=1)
326  n = n + 1
327  if (ialayer(n) == k) then
328  dbl3d(j, i, k) = dbl1d(n)
329  else
330  dbl3d(j, i, k) = dnodata
331  end if
332  end do
333  end do
334  end do
335  else
336  dbl3d = dnodata
337  nvals = export_pkg%mshape(3) * export_pkg%mshape(2)
338  dbl2d_ptr(1:export_pkg%mshape(3), 1:export_pkg%mshape(2)) => dbl1d(1:nvals)
339  dbl3d(:, :, 1) = dbl2d_ptr(:, :)
340  end if
341  !
342  call nc_export_dbl3d(this%ncid, this%dim_ids, this%var_ids, this%dis, dbl3d, &
343  nc_varname, export_pkg%mf6_input%subcomponent_name, &
344  idt%tagname, this%gridmap_name, idt%shape, &
345  idt%longname, input_attr, this%deflate, this%shuffle, &
346  this%chunk_face, export_pkg%iper, idxaux, this%nc_fname)
347 
348  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 858 of file DisNCMesh.f90.

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

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

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

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

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

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

281  class(Mesh2dDisExportType), intent(inout) :: this
282  class(ExportPackageType), pointer, intent(in) :: export_pkg
283  errmsg = 'NetCDF period export not supported for model='// &
284  trim(this%modelname)//', package='// &
285  trim(export_pkg%mf6_input%subcomponent_name)
286  call store_error(errmsg, .true.)
287  !
288  ! -- synchronize file
289  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 185 of file DisNCMesh.f90.

186  use constantsmodule, only: dnodata, dzero
187  use tdismodule, only: kper
190  class(Mesh2dDisExportType), intent(inout) :: this
191  class(ExportPackageType), pointer, intent(in) :: export_pkg
192  character(len=*), intent(in) :: ilayer_varname
193  integer(I4B), intent(in) :: ilayer
194  ! -- local
195  type(InputParamDefinitionType), pointer :: idt
196  integer(I4B), dimension(:), pointer, contiguous :: int1d
197  real(DP), dimension(:), pointer, contiguous :: dbl1d
198  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
199  integer(I4B), dimension(:), pointer, contiguous :: ialayer
200  real(DP), dimension(:), contiguous, pointer :: dbl1d_ptr
201  character(len=LINELENGTH) :: nc_varname, input_attr
202  integer(I4B) :: n, iparam, nvals
203  logical(LGP) :: ilayer_read
204  !
205  ! -- initialize
206  nullify (ialayer)
207  ilayer_read = .false.
208  !
209  ! -- set pointer to ilayer variable
210  call mem_setptr(ialayer, export_pkg%param_names(ilayer), &
211  export_pkg%mf6_input%mempath)
212  !
213  ! -- check if layer index variable was read
214  if (export_pkg%param_reads(ilayer)%invar == 1) then
215  ilayer_read = .true.
216  end if
217  !
218  ! -- export defined period input
219  do iparam = 1, export_pkg%nparam
220  !
221  ! -- check if variable was read this period
222  if (export_pkg%param_reads(iparam)%invar < 1) cycle
223  !
224  ! -- set input definition
225  idt => &
226  get_param_definition_type(export_pkg%mf6_input%param_dfns, &
227  export_pkg%mf6_input%component_type, &
228  export_pkg%mf6_input%subcomponent_type, &
229  'PERIOD', export_pkg%param_names(iparam), '')
230  !
231  ! -- set variable name and input string
232  nc_varname = trim(export_pkg%mf6_input%subcomponent_name)//'_'// &
233  trim(idt%mf6varname)
234  input_attr = this%input_attribute(export_pkg%mf6_input%subcomponent_name, &
235  idt)
236  !
237  ! -- export arrays
238  select case (idt%datatype)
239  case ('INTEGER1D')
240  call mem_setptr(int1d, idt%mf6varname, export_pkg%mf6_input%mempath)
241  call nc_export_int1d(this%ncid, this%dim_ids, this%x_dim, this%y_dim, &
242  this%var_ids, this%dis, int1d, nc_varname, &
243  export_pkg%mf6_input%subcomponent_name, &
244  idt%tagname, this%gridmap_name, idt%shape, &
245  idt%longname, input_attr, this%deflate, &
246  this%shuffle, this%chunk_face, kper, this%nc_fname)
247  case ('DOUBLE1D')
248  call mem_setptr(dbl1d, idt%mf6varname, export_pkg%mf6_input%mempath)
249  call this%export_layer_3d(export_pkg, idt, ilayer_read, ialayer, &
250  dbl1d, nc_varname, input_attr)
251  case ('DOUBLE2D')
252  call mem_setptr(dbl2d, idt%mf6varname, export_pkg%mf6_input%mempath)
253  nvals = this%dis%ncol * this%dis%nrow
254  !
255  do n = 1, size(dbl2d, dim=1) !naux
256  dbl1d_ptr(1:nvals) => dbl2d(n, :)
257  if (all(dbl1d_ptr == dzero)) then
258  else
259  call this%export_layer_3d(export_pkg, idt, ilayer_read, ialayer, &
260  dbl1d_ptr, nc_varname, input_attr, n)
261  end if
262  end do
263  case default
264  !
265  errmsg = 'EXPORT ilayaer unsupported datatype='//trim(idt%datatype)
266  call store_error(errmsg, .true.)
267  end select
268  end do
269  !
270  ! -- synchronize file
271  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
272  !
273  ! -- return
274  return
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  !
140  ! -- allocate nodereduced size 1d array
141  allocate (dbl1d(size(this%dis%nodereduced)))
142  !
143  ! -- initialize DHNOFLO for non-active cells
144  dbl1d = dhnoflo
145  !
146  ! -- update active cells
147  do n = 1, size(this%dis%nodereduced)
148  if (this%dis%nodereduced(n) > 0) then
149  dbl1d(n) = this%x(this%dis%nodereduced(n))
150  end if
151  end do
152  !
153  dbl2d(1:dis_shape(1), 1:dis_shape(2)) => dbl1d(1:nvals)
154  else
155  !
156  dbl2d(1:dis_shape(1), 1:dis_shape(2)) => this%x(1:nvals)
157  !
158  end if
159  !
160  do k = 1, this%dis%nlay
161  ! -- extend array with step data
162  call nf_verify(nf90_put_var(this%ncid, &
163  this%var_ids%dependent(k), dbl2d(:, k), &
164  start=(/1, this%stepcnt/), &
165  count=(/(this%dis%ncol * this%dis%nrow), 1/)), &
166  this%nc_fname)
167  end do
168  !
169  ! -- write to time coordinate variable
170  call nf_verify(nf90_put_var(this%ncid, this%var_ids%time, &
171  totim, start=(/this%stepcnt/)), &
172  this%nc_fname)
173  !
174  ! -- update file
175  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
176  !
177  ! --cleanup
178  if (associated(dbl1d)) deallocate (dbl1d)
179  nullify (dbl1d)
180  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: