MODFLOW 6  version 6.7.0.dev2
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 (this, export_pkg)
 netcdf export package dynamic input More...
 
subroutine export_input_array (this, pkgtype, pkgname, mempath, idt)
 netcdf export an input array More...
 
subroutine define_dim (this)
 netcdf export define dimensions More...
 
subroutine add_mesh_data (this)
 netcdf export add mesh information More...
 
subroutine nc_export_int1d (p_mem, ncid, dim_ids, x_dim, y_dim, var_ids, dis, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, iper, nc_fname)
 netcdf export 1D integer More...
 
subroutine nc_export_int2d (p_mem, ncid, dim_ids, var_ids, dis, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 2D integer More...
 
subroutine nc_export_int3d (p_mem, ncid, dim_ids, var_ids, dis, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 3D integer More...
 
subroutine nc_export_dbl1d (p_mem, ncid, dim_ids, x_dim, y_dim, var_ids, dis, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, iper, iaux, nc_fname)
 netcdf export 1D double More...
 
subroutine nc_export_dbl2d (p_mem, ncid, dim_ids, var_ids, dis, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, nc_fname)
 netcdf export 2D double More...
 
subroutine nc_export_dbl3d (p_mem, ncid, dim_ids, var_ids, dis, idt, mempath, nc_tag, pkgname, gridmap_name, deflate, shuffle, chunk_face, 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)
private

Definition at line 423 of file DisNCMesh.f90.

425  class(Mesh2dDisExportType), intent(inout) :: this
426  integer(I4B) :: cnt, maxvert, m
427  integer(I4B), dimension(:), allocatable :: verts
428  real(DP), dimension(:), allocatable :: bnds
429  integer(I4B) :: i, j
430  real(DP) :: x, y, x_transform, y_transform
431  real(DP), dimension(:), allocatable :: node_x, node_y
432  real(DP), dimension(:), allocatable :: cell_x, cell_y
433 
434  ! initialize max vertices required to define cell
435  maxvert = 4
436 
437  ! set mesh container variable value to 1
438  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh, 1), &
439  this%nc_fname)
440 
441  ! allocate temporary arrays
442  allocate (verts(maxvert))
443  allocate (bnds(maxvert))
444  allocate (node_x(((this%dis%ncol + 1) * (this%dis%nrow + 1))))
445  allocate (node_y(((this%dis%ncol + 1) * (this%dis%nrow + 1))))
446  allocate (cell_x((this%dis%ncol * this%dis%nrow)))
447  allocate (cell_y((this%dis%ncol * this%dis%nrow)))
448 
449  ! set node_x and node_y arrays
450  cnt = 0
451  node_x = nf90_fill_double
452  node_y = nf90_fill_double
453  y = sum(this%dis%delc)
454  do j = this%dis%nrow, 0, -1
455  x = 0
456  do i = this%dis%ncol, 0, -1
457  cnt = cnt + 1
458  call dis_transform_xy(x, y, &
459  this%dis%xorigin, &
460  this%dis%yorigin, &
461  this%dis%angrot, &
462  x_transform, y_transform)
463  node_x(cnt) = x_transform
464  node_y(cnt) = y_transform
465  if (i > 0) x = x + this%dis%delr(i)
466  end do
467  if (j > 0) y = y - this%dis%delc(j)
468  end do
469 
470  ! write node_x and node_y arrays to netcdf file
471  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_x, node_x), &
472  this%nc_fname)
473  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_node_y, node_y), &
474  this%nc_fname)
475 
476  ! set cell_x and cell_y arrays
477  cnt = 1
478  cell_x = nf90_fill_double
479  cell_y = nf90_fill_double
480  do j = 1, this%dis%nrow
481  y = this%dis%celly(j)
482  do i = 1, this%dis%ncol
483  x = this%dis%cellx(i)
484  call dis_transform_xy(x, y, &
485  this%dis%xorigin, &
486  this%dis%yorigin, &
487  this%dis%angrot, &
488  x_transform, y_transform)
489  cell_x(cnt) = x_transform
490  cell_y(cnt) = y_transform
491  cnt = cnt + 1
492  end do
493  end do
494 
495  ! write face_x and face_y arrays to netcdf file
496  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_x, cell_x), &
497  this%nc_fname)
498  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_y, cell_y), &
499  this%nc_fname)
500 
501  ! set face nodes array
502  cnt = 0
503  do i = 1, this%dis%nrow
504  do j = 1, this%dis%ncol
505  cnt = cnt + 1
506  verts = nf90_fill_int
507  verts(1) = cnt + this%dis%ncol + i
508  verts(2) = cnt + this%dis%ncol + i + 1
509  if (i > 1) then
510  verts(3) = cnt + i
511  verts(4) = cnt + i - 1
512  else
513  verts(3) = cnt + 1
514  verts(4) = cnt
515  end if
516 
517  ! write face nodes array to netcdf file
518  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_nodes, &
519  verts, start=(/1, cnt/), &
520  count=(/maxvert, 1/)), &
521  this%nc_fname)
522 
523  ! set face y bounds array
524  bnds = nf90_fill_double
525  do m = 1, size(bnds)
526  if (verts(m) /= nf90_fill_int) then
527  bnds(m) = node_y(verts(m))
528  end if
529  ! write face y bounds array to netcdf file
530  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_ybnds, &
531  bnds, start=(/1, cnt/), &
532  count=(/maxvert, 1/)), &
533  this%nc_fname)
534  end do
535 
536  ! set face x bounds array
537  bnds = nf90_fill_double
538  do m = 1, size(bnds)
539  if (verts(m) /= nf90_fill_int) then
540  bnds(m) = node_x(verts(m))
541  end if
542  ! write face x bounds array to netcdf file
543  call nf_verify(nf90_put_var(this%ncid, this%var_ids%mesh_face_xbnds, &
544  bnds, start=(/1, cnt/), &
545  count=(/maxvert, 1/)), &
546  this%nc_fname)
547  end do
548  end do
549  end do
550 
551  ! cleanup
552  deallocate (bnds)
553  deallocate (verts)
554  deallocate (node_x)
555  deallocate (node_y)
556  deallocate (cell_x)
557  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 381 of file DisNCMesh.f90.

382  class(Mesh2dDisExportType), intent(inout) :: this
383 
384  ! time
385  call nf_verify(nf90_def_dim(this%ncid, 'time', this%totnstp, &
386  this%dim_ids%time), this%nc_fname)
387  call nf_verify(nf90_def_var(this%ncid, 'time', nf90_double, &
388  this%dim_ids%time, this%var_ids%time), &
389  this%nc_fname)
390  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'calendar', &
391  'standard'), this%nc_fname)
392  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'units', &
393  this%datetime), this%nc_fname)
394  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'axis', 'T'), &
395  this%nc_fname)
396  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'standard_name', &
397  'time'), this%nc_fname)
398  call nf_verify(nf90_put_att(this%ncid, this%var_ids%time, 'long_name', &
399  'time'), this%nc_fname)
400 
401  ! mesh
402  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_node', &
403  ((this%dis%ncol + 1) * (this%dis%nrow + 1)), &
404  this%dim_ids%nmesh_node), this%nc_fname)
405  call nf_verify(nf90_def_dim(this%ncid, 'nmesh_face', &
406  (this%dis%ncol * this%dis%nrow), &
407  this%dim_ids%nmesh_face), this%nc_fname)
408  call nf_verify(nf90_def_dim(this%ncid, 'max_nmesh_face_nodes', 4, &
409  this%dim_ids%max_nmesh_face_nodes), &
410  this%nc_fname)
411 
412  ! x, y, nlay
413  call nf_verify(nf90_def_dim(this%ncid, 'nlay', this%dis%nlay, &
414  this%dim_ids%nlay), this%nc_fname)
415  call nf_verify(nf90_def_dim(this%ncid, 'x', this%dis%ncol, &
416  this%x_dim), this%nc_fname)
417  call nf_verify(nf90_def_dim(this%ncid, 'y', this%dis%nrow, &
418  this%y_dim), this%nc_fname)
Here is the call graph for this function:

◆ df()

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

Definition at line 86 of file DisNCMesh.f90.

87  use constantsmodule, only: mvalidate
88  use simvariablesmodule, only: isim_mode
89  class(Mesh2dDisExportType), intent(inout) :: this
90  ! put root group file scope attributes
91  call this%add_global_att()
92  ! define root group dimensions and coordinate variables
93  call this%define_dim()
94  ! define mesh variables
95  call this%create_mesh()
96  if (isim_mode /= mvalidate) then
97  ! define the dependent variable
98  call this%define_dependent()
99  end if
100  ! define period input arrays
101  call this%df_export()
102  ! exit define mode
103  call nf_verify(nf90_enddef(this%ncid), this%nc_fname)
104  ! create mesh
105  call this%add_mesh_data()
106  ! define and set package input griddata
107  call this%add_pkg_data()
108  ! define and set gridmap variable
109  call this%define_gridmap()
110  ! synchronize file
111  call nf_verify(nf90_sync(this%ncid), 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:

◆ dis_export_destroy()

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

Definition at line 76 of file DisNCMesh.f90.

77  class(Mesh2dDisExportType), intent(inout) :: this
78  deallocate (this%var_ids%dependent)
79  ! destroy base class
80  call this%mesh_destroy()
81  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 50 of file DisNCMesh.f90.

53  class(Mesh2dDisExportType), intent(inout) :: this
54  character(len=*), intent(in) :: modelname
55  character(len=*), intent(in) :: modeltype
56  character(len=*), intent(in) :: modelfname
57  character(len=*), intent(in) :: nc_fname
58  integer(I4B), intent(in) :: disenum
59  integer(I4B), intent(in) :: nctype
60  integer(I4B), intent(in) :: iout
61 
62  ! set nlay
63  this%nlay = this%dis%nlay
64 
65  ! allocate var_id arrays
66  allocate (this%var_ids%dependent(this%nlay))
67  allocate (this%var_ids%export(this%nlay))
68 
69  ! initialize base class
70  call this%mesh_init(modelname, modeltype, modelfname, nc_fname, disenum, &
71  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 315 of file DisNCMesh.f90.

316  class(Mesh2dDisExportType), intent(inout) :: this
317  character(len=*), intent(in) :: pkgtype
318  character(len=*), intent(in) :: pkgname
319  character(len=*), intent(in) :: mempath
320  type(InputParamDefinitionType), pointer, intent(in) :: idt
321  integer(I4B), dimension(:), pointer, contiguous :: int1d
322  integer(I4B), dimension(:, :), pointer, contiguous :: int2d
323  integer(I4B), dimension(:, :, :), pointer, contiguous :: int3d
324  real(DP), dimension(:), pointer, contiguous :: dbl1d
325  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
326  real(DP), dimension(:, :, :), pointer, contiguous :: dbl3d
327  character(len=LINELENGTH) :: nc_tag
328  integer(I4B) :: iper, iaux
329 
330  iper = 0
331  iaux = 0
332 
333  ! set package input tag
334  nc_tag = this%input_attribute(pkgname, idt)
335 
336  select case (idt%datatype)
337  case ('INTEGER1D')
338  call mem_setptr(int1d, idt%mf6varname, mempath)
339  call nc_export_int1d(int1d, this%ncid, this%dim_ids, this%x_dim, &
340  this%y_dim, this%var_ids, this%dis, idt, mempath, &
341  nc_tag, pkgname, this%gridmap_name, this%deflate, &
342  this%shuffle, this%chunk_face, iper, this%nc_fname)
343  case ('INTEGER2D')
344  call mem_setptr(int2d, idt%mf6varname, mempath)
345  call nc_export_int2d(int2d, this%ncid, this%dim_ids, this%var_ids, &
346  this%dis, idt, mempath, nc_tag, pkgname, &
347  this%gridmap_name, this%deflate, this%shuffle, &
348  this%chunk_face, this%nc_fname)
349  case ('INTEGER3D')
350  call mem_setptr(int3d, idt%mf6varname, mempath)
351  call nc_export_int3d(int3d, this%ncid, this%dim_ids, this%var_ids, &
352  this%dis, idt, mempath, nc_tag, pkgname, &
353  this%gridmap_name, this%deflate, this%shuffle, &
354  this%chunk_face, this%nc_fname)
355  case ('DOUBLE1D')
356  call mem_setptr(dbl1d, idt%mf6varname, mempath)
357  call nc_export_dbl1d(dbl1d, this%ncid, this%dim_ids, this%x_dim, &
358  this%y_dim, this%var_ids, this%dis, idt, mempath, &
359  nc_tag, pkgname, this%gridmap_name, this%deflate, &
360  this%shuffle, this%chunk_face, iper, iaux, &
361  this%nc_fname)
362  case ('DOUBLE2D')
363  call mem_setptr(dbl2d, idt%mf6varname, mempath)
364  call nc_export_dbl2d(dbl2d, this%ncid, this%dim_ids, this%var_ids, &
365  this%dis, idt, mempath, nc_tag, pkgname, &
366  this%gridmap_name, this%deflate, this%shuffle, &
367  this%chunk_face, this%nc_fname)
368  case ('DOUBLE3D')
369  call mem_setptr(dbl3d, idt%mf6varname, mempath)
370  call nc_export_dbl3d(dbl3d, this%ncid, this%dim_ids, this%var_ids, &
371  this%dis, idt, mempath, nc_tag, pkgname, &
372  this%gridmap_name, this%deflate, this%shuffle, &
373  this%chunk_face, this%nc_fname)
374  case default
375  ! no-op, no other datatypes exported
376  end select
Here is the call graph for this function:

◆ nc_export_dbl1d()

subroutine meshdismodelmodule::nc_export_dbl1d ( real(dp), dimension(:), intent(in), pointer, contiguous  p_mem,
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(in)  var_ids,
type(distype), intent(in), pointer  dis,
type(inputparamdefinitiontype), pointer  idt,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  gridmap_name,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
integer(i4b), intent(in)  iper,
integer(i4b), intent(in)  iaux,
character(len=*), intent(in)  nc_fname 
)
private

Definition at line 828 of file DisNCMesh.f90.

831  use netcdfcommonmodule, only: ixstp
832  real(DP), dimension(:), pointer, contiguous, intent(in) :: p_mem
833  integer(I4B), intent(in) :: ncid
834  type(MeshNCDimIdType), intent(inout) :: dim_ids
835  integer(I4B), intent(in) :: x_dim
836  integer(I4B), intent(in) :: y_dim
837  type(MeshNCVarIdType), intent(in) :: var_ids
838  type(DisType), pointer, intent(in) :: dis
839  type(InputParamDefinitionType), pointer :: idt
840  character(len=*), intent(in) :: mempath
841  character(len=*), intent(in) :: nc_tag
842  character(len=*), intent(in) :: pkgname
843  character(len=*), intent(in) :: gridmap_name
844  integer(I4B), intent(in) :: deflate
845  integer(I4B), intent(in) :: shuffle
846  integer(I4B), intent(in) :: chunk_face
847  integer(I4B), intent(in) :: iper
848  integer(I4B), intent(in) :: iaux
849  character(len=*), intent(in) :: nc_fname
850  real(DP), dimension(:, :, :), pointer, contiguous :: dbl3d
851  real(DP), dimension(:), pointer, contiguous :: dbl1d
852  integer(I4B) :: axis_dim, nvals, k, istp
853  integer(NF90_INT), dimension(:), allocatable :: var_id
854  character(len=LINELENGTH) :: longname, varname
855 
856  if (idt%shape == 'NROW' .or. &
857  idt%shape == 'NCOL' .or. &
858  idt%shape == 'NCPL' .or. &
859  idt%shape == 'NAUX NCPL') then
860 
861  if (iper == 0) then
862 
863  select case (idt%shape)
864  case ('NROW')
865  axis_dim = y_dim
866  case ('NCOL')
867  axis_dim = x_dim
868  case ('NCPL', 'NAUX NCPL')
869  axis_dim = dim_ids%nmesh_face
870  end select
871 
872  ! set names
873  varname = export_varname(pkgname, idt%tagname, mempath, iaux=iaux)
874  longname = export_longname(idt%longname, pkgname, idt%tagname, &
875  mempath, iaux=iaux)
876 
877  allocate (var_id(1))
878 
879  ! reenter define mode and create variable
880  call nf_verify(nf90_redef(ncid), nc_fname)
881  call nf_verify(nf90_def_var(ncid, varname, nf90_double, &
882  (/axis_dim/), var_id(1)), &
883  nc_fname)
884 
885  ! NROW/NCOL shapes use default chunking
886  call ncvar_deflate(ncid, var_id(1), deflate, shuffle, nc_fname)
887 
888  ! put attr
889  call nf_verify(nf90_put_att(ncid, var_id(1), '_FillValue', &
890  (/nf90_fill_double/)), nc_fname)
891  call nf_verify(nf90_put_att(ncid, var_id(1), 'long_name', &
892  longname), nc_fname)
893 
894  ! add mf6 attr
895  call ncvar_mf6attr(ncid, var_id(1), 0, iaux, nc_tag, nc_fname)
896 
897  ! exit define mode and write data
898  call nf_verify(nf90_enddef(ncid), nc_fname)
899  call nf_verify(nf90_put_var(ncid, var_id(1), p_mem), &
900  nc_fname)
901  else
902  istp = ixstp()
903  nvals = dis%nrow * dis%ncol
904  call nf_verify(nf90_put_var(ncid, &
905  var_ids%export(1), p_mem, &
906  start=(/1, istp/), &
907  count=(/nvals, 1/)), nc_fname)
908  end if
909 
910  else
911  ! reshape input
912  dbl3d(1:dis%ncol, 1:dis%nrow, 1:dis%nlay) => p_mem(1:dis%nodesuser)
913 
914  ! set nvals as ncpl
915  nvals = dis%nrow * dis%ncol
916 
917  if (iper == 0) then
918  ! not a timeseries, create variables and write griddata
919 
920  ! allocate local variable id storage
921  allocate (var_id(dis%nlay))
922 
923  ! reenter define mode and create layer variables
924  call nf_verify(nf90_redef(ncid), nc_fname)
925  do k = 1, dis%nlay
926  ! set names
927  varname = export_varname(pkgname, idt%tagname, mempath, layer=k, &
928  iaux=iaux)
929  longname = export_longname(idt%longname, pkgname, idt%tagname, &
930  mempath, layer=k, iaux=iaux)
931 
932  ! create layer variable
933  call nf_verify(nf90_def_var(ncid, varname, nf90_double, &
934  (/dim_ids%nmesh_face/), var_id(k)), &
935  nc_fname)
936 
937  ! apply chunking parameters
938  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
939  ! deflate and shuffle
940  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
941 
942  ! put attr
943  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
944  (/nf90_fill_double/)), nc_fname)
945  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
946  longname), nc_fname)
947 
948  ! add grid mapping and mf6 attr
949  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
950  call ncvar_mf6attr(ncid, var_id(k), k, iaux, nc_tag, nc_fname)
951  end do
952 
953  ! exit define mode
954  call nf_verify(nf90_enddef(ncid), nc_fname)
955 
956  ! write layer data
957  do k = 1, dis%nlay
958  dbl1d(1:nvals) => dbl3d(:, :, k)
959  call nf_verify(nf90_put_var(ncid, var_id(k), dbl1d), nc_fname)
960  end do
961 
962  ! cleanup
963  deallocate (var_id)
964  else
965  ! timeseries, add period data
966  istp = ixstp()
967  do k = 1, dis%nlay
968  dbl1d(1:nvals) => dbl3d(:, :, k)
969  call nf_verify(nf90_put_var(ncid, &
970  var_ids%export(k), dbl1d, &
971  start=(/1, istp/), &
972  count=(/nvals, 1/)), nc_fname)
973  end do
974  end if
975  end if
This module contains the NetCDFCommonModule.
Definition: NetCDFCommon.f90:6
integer(i4b) function, public ixstp()
step index for timeseries data
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nc_export_dbl2d()

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

Definition at line 980 of file DisNCMesh.f90.

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

Definition at line 1035 of file DisNCMesh.f90.

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

Definition at line 562 of file DisNCMesh.f90.

565  use netcdfcommonmodule, only: ixstp
566  integer(I4B), dimension(:), pointer, contiguous, intent(in) :: p_mem
567  integer(I4B), intent(in) :: ncid
568  type(MeshNCDimIdType), intent(inout) :: dim_ids
569  integer(I4B), intent(in) :: x_dim
570  integer(I4B), intent(in) :: y_dim
571  type(MeshNCVarIdType), intent(inout) :: var_ids
572  type(DisType), pointer, intent(in) :: dis
573  type(InputParamDefinitionType), pointer :: idt
574  character(len=*), intent(in) :: mempath
575  character(len=*), intent(in) :: nc_tag
576  character(len=*), intent(in) :: pkgname
577  character(len=*), intent(in) :: gridmap_name
578  integer(I4B), intent(in) :: deflate
579  integer(I4B), intent(in) :: shuffle
580  integer(I4B), intent(in) :: chunk_face
581  integer(I4B), intent(in) :: iper
582  character(len=*), intent(in) :: nc_fname
583  integer(I4B), dimension(:, :, :), pointer, contiguous :: int3d
584  integer(I4B), dimension(:), pointer, contiguous :: int1d
585  integer(I4B) :: axis_dim, nvals, k, istp
586  integer(I4B), dimension(:), allocatable :: var_id
587  character(len=LINELENGTH) :: longname, varname
588 
589  if (idt%shape == 'NROW' .or. &
590  idt%shape == 'NCOL' .or. &
591  idt%shape == 'NCPL' .or. &
592  idt%shape == 'NAUX NCPL') then
593 
594  if (iper == 0) then
595 
596  select case (idt%shape)
597  case ('NROW')
598  axis_dim = y_dim
599  case ('NCOL')
600  axis_dim = x_dim
601  case ('NCPL', 'NAUX NCPL')
602  axis_dim = dim_ids%nmesh_face
603  end select
604 
605  ! set names
606  varname = export_varname(pkgname, idt%tagname, mempath)
607  longname = export_longname(idt%longname, pkgname, idt%tagname, mempath)
608 
609  allocate (var_id(1))
610 
611  ! reenter define mode and create variable
612  call nf_verify(nf90_redef(ncid), nc_fname)
613  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
614  (/axis_dim/), var_id(1)), &
615  nc_fname)
616 
617  ! NROW/NCOL shapes use default chunking
618  call ncvar_deflate(ncid, var_id(1), deflate, shuffle, nc_fname)
619 
620  ! put attr
621  call nf_verify(nf90_put_att(ncid, var_id(1), '_FillValue', &
622  (/nf90_fill_int/)), nc_fname)
623  call nf_verify(nf90_put_att(ncid, var_id(1), 'long_name', &
624  longname), nc_fname)
625 
626  ! add mf6 attr
627  call ncvar_mf6attr(ncid, var_id(1), 0, 0, nc_tag, nc_fname)
628 
629  ! exit define mode and write data
630  call nf_verify(nf90_enddef(ncid), nc_fname)
631  call nf_verify(nf90_put_var(ncid, var_id(1), p_mem), &
632  nc_fname)
633  else
634  istp = ixstp()
635  nvals = dis%nrow * dis%ncol
636  call nf_verify(nf90_put_var(ncid, &
637  var_ids%export(1), p_mem, &
638  start=(/1, istp/), &
639  count=(/nvals, 1/)), nc_fname)
640  end if
641 
642  else
643  ! reshape input
644  int3d(1:dis%ncol, 1:dis%nrow, 1:dis%nlay) => p_mem(1:dis%nodesuser)
645 
646  ! set nvals as ncpl
647  nvals = dis%nrow * dis%ncol
648 
649  if (iper == 0) then
650  ! not a timeseries, create variables and write griddata
651  allocate (var_id(dis%nlay))
652 
653  ! reenter define mode and create variable
654  call nf_verify(nf90_redef(ncid), nc_fname)
655  do k = 1, dis%nlay
656  ! set names
657  varname = export_varname(pkgname, idt%tagname, mempath, &
658  layer=k)
659  longname = export_longname(idt%longname, pkgname, idt%tagname, &
660  mempath, layer=k)
661 
662  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
663  (/dim_ids%nmesh_face/), var_id(k)), &
664  nc_fname)
665 
666  ! apply chunking parameters
667  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
668  ! deflate and shuffle
669  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
670 
671  ! put attr
672  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
673  (/nf90_fill_int/)), nc_fname)
674  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
675  longname), nc_fname)
676 
677  ! add grid mapping and mf6 attr
678  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
679  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
680  end do
681 
682  ! exit define mode and write data
683  call nf_verify(nf90_enddef(ncid), nc_fname)
684  do k = 1, dis%nlay
685  int1d(1:nvals) => int3d(:, :, k)
686  call nf_verify(nf90_put_var(ncid, var_id(k), int1d), nc_fname)
687  end do
688 
689  ! cleanup
690  deallocate (var_id)
691  else
692  ! timeseries, add period data
693  istp = ixstp()
694  do k = 1, dis%nlay
695  int1d(1:nvals) => int3d(:, :, k)
696  call nf_verify(nf90_put_var(ncid, &
697  var_ids%export(k), int1d, &
698  start=(/1, istp/), &
699  count=(/nvals, 1/)), nc_fname)
700  end do
701  end if
702  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), dimension(:, :), intent(in), pointer, contiguous  p_mem,
integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
type(inputparamdefinitiontype), pointer  idt,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  gridmap_name,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)

Definition at line 707 of file DisNCMesh.f90.

710  integer(I4B), dimension(:, :), pointer, contiguous, intent(in) :: p_mem
711  integer(I4B), intent(in) :: ncid
712  type(MeshNCDimIdType), intent(inout) :: dim_ids
713  type(MeshNCVarIdType), intent(inout) :: var_ids
714  type(DisType), pointer, intent(in) :: dis
715  type(InputParamDefinitionType), pointer :: idt
716  character(len=*), intent(in) :: mempath
717  character(len=*), intent(in) :: nc_tag
718  character(len=*), intent(in) :: pkgname
719  character(len=*), intent(in) :: gridmap_name
720  integer(I4B), intent(in) :: deflate
721  integer(I4B), intent(in) :: shuffle
722  integer(I4B), intent(in) :: chunk_face
723  character(len=*), intent(in) :: nc_fname
724  integer(I4B) :: var_id, nvals
725  integer(I4B), dimension(:), pointer, contiguous :: int1d
726  character(len=LINELENGTH) :: longname, varname
727 
728  ! set names
729  varname = export_varname(pkgname, idt%tagname, mempath)
730  longname = export_longname(idt%longname, pkgname, idt%tagname, mempath)
731 
732  ! reenter define mode and create variable
733  call nf_verify(nf90_redef(ncid), nc_fname)
734  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
735  (/dim_ids%nmesh_face/), var_id), &
736  nc_fname)
737 
738  ! apply chunking parameters
739  call ncvar_chunk(ncid, var_id, chunk_face, nc_fname)
740  ! deflate and shuffle
741  call ncvar_deflate(ncid, var_id, deflate, shuffle, nc_fname)
742 
743  ! put attr
744  call nf_verify(nf90_put_att(ncid, var_id, '_FillValue', &
745  (/nf90_fill_int/)), nc_fname)
746  call nf_verify(nf90_put_att(ncid, var_id, 'long_name', &
747  longname), nc_fname)
748 
749  ! add grid mapping and mf6 attr
750  call ncvar_gridmap(ncid, var_id, gridmap_name, nc_fname)
751  call ncvar_mf6attr(ncid, var_id, 0, 0, nc_tag, nc_fname)
752 
753  ! exit define mode and write data
754  call nf_verify(nf90_enddef(ncid), nc_fname)
755  nvals = dis%nrow * dis%ncol
756  int1d(1:nvals) => p_mem
757  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), dimension(:, :, :), intent(in), pointer, contiguous  p_mem,
integer(i4b), intent(in)  ncid,
type(meshncdimidtype), intent(inout)  dim_ids,
type(meshncvaridtype), intent(inout)  var_ids,
type(distype), intent(in), pointer  dis,
type(inputparamdefinitiontype), pointer  idt,
character(len=*), intent(in)  mempath,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  gridmap_name,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)
private

Definition at line 762 of file DisNCMesh.f90.

765  integer(I4B), dimension(:, :, :), pointer, contiguous, intent(in) :: p_mem
766  integer(I4B), intent(in) :: ncid
767  type(MeshNCDimIdType), intent(inout) :: dim_ids
768  type(MeshNCVarIdType), intent(inout) :: var_ids
769  type(DisType), pointer, intent(in) :: dis
770  type(InputParamDefinitionType), pointer :: idt
771  character(len=*), intent(in) :: mempath
772  character(len=*), intent(in) :: nc_tag
773  character(len=*), intent(in) :: pkgname
774  character(len=*), intent(in) :: gridmap_name
775  integer(I4B), intent(in) :: deflate
776  integer(I4B), intent(in) :: shuffle
777  integer(I4B), intent(in) :: chunk_face
778  character(len=*), intent(in) :: nc_fname
779  integer(I4B), dimension(:), allocatable :: var_id
780  integer(I4B), dimension(:), pointer, contiguous :: int1d
781  character(len=LINELENGTH) :: longname, varname
782  integer(I4B) :: k, nvals
783 
784  allocate (var_id(dis%nlay))
785 
786  ! reenter define mode and create variable
787  call nf_verify(nf90_redef(ncid), nc_fname)
788  do k = 1, dis%nlay
789  ! set names
790  varname = export_varname(pkgname, idt%tagname, mempath, layer=k)
791  longname = export_longname(idt%longname, pkgname, idt%tagname, &
792  mempath, layer=k)
793 
794  call nf_verify(nf90_def_var(ncid, varname, nf90_int, &
795  (/dim_ids%nmesh_face/), var_id(k)), &
796  nc_fname)
797 
798  ! apply chunking parameters
799  call ncvar_chunk(ncid, var_id(k), chunk_face, nc_fname)
800  ! deflate and shuffle
801  call ncvar_deflate(ncid, var_id(k), deflate, shuffle, nc_fname)
802 
803  ! put attr
804  call nf_verify(nf90_put_att(ncid, var_id(k), '_FillValue', &
805  (/nf90_fill_int/)), nc_fname)
806  call nf_verify(nf90_put_att(ncid, var_id(k), 'long_name', &
807  longname), nc_fname)
808 
809  ! add grid mapping and mf6 attr
810  call ncvar_gridmap(ncid, var_id(k), gridmap_name, nc_fname)
811  call ncvar_mf6attr(ncid, var_id(k), k, 0, nc_tag, nc_fname)
812  end do
813 
814  ! exit define mode and write data
815  call nf_verify(nf90_enddef(ncid), nc_fname)
816  nvals = dis%nrow * dis%ncol
817  do k = 1, dis%nlay
818  int1d(1:nvals) => p_mem(:, :, k)
819  call nf_verify(nf90_put_var(ncid, var_id(k), int1d), nc_fname)
820  end do
821 
822  ! cleanup
823  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 183 of file DisNCMesh.f90.

184  use tdismodule, only: kper
187  class(Mesh2dDisExportType), intent(inout) :: this
188  class(ExportPackageType), pointer, intent(in) :: export_pkg
189  type(InputParamDefinitionType), pointer :: idt
190  integer(I4B), dimension(:), pointer, contiguous :: int1d
191  real(DP), dimension(:), pointer, contiguous :: dbl1d, nodes
192  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
193  character(len=LINELENGTH) :: nc_tag
194  integer(I4B) :: iaux, iparam, nvals
195  integer(I4B) :: k, n
196  integer(I4B), pointer :: nbound
197 
198  ! initialize
199  iaux = 0
200 
201  ! export defined period input
202  do iparam = 1, export_pkg%nparam
203  ! check if variable was read this period
204  if (export_pkg%param_reads(iparam)%invar < 1) cycle
205 
206  ! set input definition
207  idt => &
208  get_param_definition_type(export_pkg%mf6_input%param_dfns, &
209  export_pkg%mf6_input%component_type, &
210  export_pkg%mf6_input%subcomponent_type, &
211  'PERIOD', export_pkg%param_names(iparam), '')
212 
213  ! set variable input tag
214  nc_tag = this%input_attribute(export_pkg%mf6_input%subcomponent_name, &
215  idt)
216 
217  ! export arrays
218  select case (idt%datatype)
219  case ('INTEGER1D')
220  call mem_setptr(int1d, idt%mf6varname, export_pkg%mf6_input%mempath)
221  this%var_ids%export(1) = export_pkg%varids_param(iparam, 1)
222  call nc_export_int1d(int1d, this%ncid, this%dim_ids, this%x_dim, &
223  this%y_dim, this%var_ids, this%dis, idt, &
224  export_pkg%mf6_input%mempath, nc_tag, &
225  export_pkg%mf6_input%subcomponent_name, &
226  this%gridmap_name, this%deflate, this%shuffle, &
227  this%chunk_face, kper, this%nc_fname)
228  case ('DOUBLE1D')
229  call mem_setptr(dbl1d, idt%mf6varname, export_pkg%mf6_input%mempath)
230  select case (idt%shape)
231  case ('NCPL')
232  this%var_ids%export(1) = export_pkg%varids_param(iparam, 1)
233  call nc_export_dbl1d(dbl1d, this%ncid, this%dim_ids, this%x_dim, &
234  this%y_dim, this%var_ids, this%dis, idt, &
235  export_pkg%mf6_input%mempath, nc_tag, &
236  export_pkg%mf6_input%subcomponent_name, &
237  this%gridmap_name, this%deflate, this%shuffle, &
238  this%chunk_face, kper, iaux, this%nc_fname)
239  case ('NODES')
240  nvals = this%dis%nodesuser
241  allocate (nodes(nvals))
242  nodes = dnodata
243  do k = 1, this%dis%nlay
244  this%var_ids%export(k) = export_pkg%varids_param(iparam, k)
245  end do
246  call mem_setptr(dbl1d, idt%mf6varname, export_pkg%mf6_input%mempath)
247  call mem_setptr(int1d, 'NODEULIST', export_pkg%mf6_input%mempath)
248  call mem_setptr(nbound, 'NBOUND', export_pkg%mf6_input%mempath)
249  do n = 1, nbound
250  nodes(int1d(n)) = dbl1d(n)
251  end do
252  call nc_export_dbl1d(nodes, this%ncid, this%dim_ids, this%x_dim, &
253  this%y_dim, this%var_ids, this%dis, idt, &
254  export_pkg%mf6_input%mempath, nc_tag, &
255  export_pkg%mf6_input%subcomponent_name, &
256  this%gridmap_name, this%deflate, this%shuffle, &
257  this%chunk_face, kper, iaux, this%nc_fname)
258  deallocate (nodes)
259  case default
260  end select
261  case ('DOUBLE2D')
262  call mem_setptr(dbl2d, idt%mf6varname, export_pkg%mf6_input%mempath)
263  select case (idt%shape)
264  case ('NAUX NCPL')
265  nvals = this%dis%nrow * this%dis%ncol
266  allocate (nodes(nvals))
267  do iaux = 1, size(dbl2d, dim=1) !naux
268  this%var_ids%export(1) = export_pkg%varids_aux(iaux, 1)
269  do n = 1, nvals
270  nodes(n) = dbl2d(iaux, n)
271  end do
272  call nc_export_dbl1d(nodes, this%ncid, this%dim_ids, this%x_dim, &
273  this%y_dim, this%var_ids, this%dis, idt, &
274  export_pkg%mf6_input%mempath, nc_tag, &
275  export_pkg%mf6_input%subcomponent_name, &
276  this%gridmap_name, this%deflate, this%shuffle, &
277  this%chunk_face, kper, iaux, this%nc_fname)
278  end do
279  deallocate (nodes)
280  case ('NAUX NODES')
281  nvals = this%dis%nodesuser
282  allocate (nodes(nvals))
283  call mem_setptr(int1d, 'NODEULIST', export_pkg%mf6_input%mempath)
284  call mem_setptr(nbound, 'NBOUND', export_pkg%mf6_input%mempath)
285  do iaux = 1, size(dbl2d, dim=1) ! naux
286  nodes = dnodata
287  do k = 1, this%dis%nlay
288  this%var_ids%export(k) = export_pkg%varids_aux(iaux, k)
289  end do
290  do n = 1, nbound
291  nodes(int1d(n)) = dbl2d(iaux, n)
292  end do
293  call nc_export_dbl1d(nodes, this%ncid, this%dim_ids, this%x_dim, &
294  this%y_dim, this%var_ids, this%dis, idt, &
295  export_pkg%mf6_input%mempath, nc_tag, &
296  export_pkg%mf6_input%subcomponent_name, &
297  this%gridmap_name, this%deflate, this%shuffle, &
298  this%chunk_face, kper, iaux, this%nc_fname)
299 
300  end do
301  deallocate (nodes)
302  case default
303  end select
304  case default
305  ! no-op, no other datatypes exported
306  end select
307  end do
308 
309  ! synchronize file
310  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
This module contains the DefinitionSelectModule.
type(inputparamdefinitiontype) function, pointer, public get_param_definition_type(input_definition_types, component_type, subcomponent_type, blockname, tagname, filename)
Return parameter definition.
This module contains the NCModelExportModule.
Definition: NCModel.f90:8
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 116 of file DisNCMesh.f90.

117  use constantsmodule, only: dhnoflo
118  use tdismodule, only: totim
119  use netcdfcommonmodule, only: ixstp
120  class(Mesh2dDisExportType), intent(inout) :: this
121  real(DP), dimension(:), pointer, contiguous :: dbl1d
122  integer(I4B) :: n, k, nvals, istp
123  integer(I4B), dimension(2) :: dis_shape
124  real(DP), dimension(:, :), pointer, contiguous :: dbl2d
125 
126  ! initialize
127  nullify (dbl1d)
128  nullify (dbl2d)
129 
130  ! set global step index
131  istp = ixstp()
132 
133  dis_shape(1) = this%dis%ncol * this%dis%nrow
134  dis_shape(2) = this%dis%nlay
135 
136  nvals = product(dis_shape)
137 
138  ! add data to dependent variable
139  if (size(this%dis%nodeuser) < &
140  size(this%dis%nodereduced)) then
141  ! allocate nodereduced size 1d array
142  allocate (dbl1d(size(this%dis%nodereduced)))
143 
144  ! initialize DHNOFLO for non-active cells
145  dbl1d = dhnoflo
146 
147  ! update active cells
148  do n = 1, size(this%dis%nodereduced)
149  if (this%dis%nodereduced(n) > 0) then
150  dbl1d(n) = this%x(this%dis%nodereduced(n))
151  end if
152  end do
153 
154  dbl2d(1:dis_shape(1), 1:dis_shape(2)) => dbl1d(1:nvals)
155  else
156  dbl2d(1:dis_shape(1), 1:dis_shape(2)) => this%x(1:nvals)
157  end if
158 
159  do k = 1, this%dis%nlay
160  ! extend array with step data
161  call nf_verify(nf90_put_var(this%ncid, &
162  this%var_ids%dependent(k), dbl2d(:, k), &
163  start=(/1, istp/), &
164  count=(/(this%dis%ncol * this%dis%nrow), 1/)), &
165  this%nc_fname)
166  end do
167 
168  ! write to time coordinate variable
169  call nf_verify(nf90_put_var(this%ncid, this%var_ids%time, &
170  totim, start=(/istp/)), &
171  this%nc_fname)
172  ! update file
173  call nf_verify(nf90_sync(this%ncid), this%nc_fname)
174 
175  ! cleanup
176  if (associated(dbl1d)) deallocate (dbl1d)
177  nullify (dbl1d)
178  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: