MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
meshmodelmodule Module Reference

This module contains the MeshModelModule. More...

Data Types

type  meshncdimidtype
 type for storing model export dimension ids More...
 
type  meshncvaridtype
 type for storing model export variable ids More...
 
type  meshmodeltype
 base ugrid netcdf export type More...
 
interface  nc_array_export_if
 abstract interfaces for derived ugrid netcd export types More...
 
type  mesh2dmodeltype
 

Functions/Subroutines

subroutine mesh_init (this, modelname, modeltype, modelfname, nc_fname, disenum, nctype, lenuni, iout)
 initialize More...
 
subroutine mesh_destroy (this)
 initialize More...
 
subroutine df_export (this)
 define timeseries input variables More...
 
subroutine export_df (this, export_pkg)
 define export package More...
 
subroutine create_timeseries (this, idt, iparam, iaux, layer, export_pkg)
 create timeseries export variable More...
 
subroutine add_global_att (this)
 create file (group) attributes More...
 
subroutine export_input_arrays (this, pkgtype, pkgname, mempath, param_dfns)
 write package gridded input data More...
 
subroutine add_pkg_data (this)
 determine packages to write gridded input More...
 
subroutine define_dependent (this)
 create the model layer dependent variables More...
 
subroutine define_gridmap (this)
 create the file grid mapping container variable More...
 
subroutine create_mesh (this)
 create the file mesh container variable More...
 
subroutine, public ncvar_chunk (ncid, varid, chunk_face, nc_fname)
 define variable chunking More...
 
subroutine, public ncvar_deflate (ncid, varid, deflate, shuffle, nc_fname)
 define variable compression More...
 
subroutine, public ncvar_gridmap (ncid, varid, gridmap_name, nc_fname)
 put variable gridmap attributes More...
 
subroutine, public ncvar_mf6attr (ncid, varid, layer, iaux, nc_tag, nc_fname)
 put variable internal attributes More...
 

Detailed Description

This module defines a base class for UGRID based (mesh) model netcdf exports. It is dependent on external netcdf libraries.

Function/Subroutine Documentation

◆ add_global_att()

subroutine meshmodelmodule::add_global_att ( class(meshmodeltype), intent(inout)  this)

Definition at line 315 of file MeshNCModel.f90.

316  class(MeshModelType), intent(inout) :: this
317  ! file scoped title
318  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'title', &
319  this%annotation%title), this%nc_fname)
320  ! source (MODFLOW 6)
321  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'source', &
322  this%annotation%source), this%nc_fname)
323  ! grid type (MODFLOW 6)
324  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'modflow_grid', &
325  this%annotation%grid), this%nc_fname)
326  ! mesh type (MODFLOW 6)
327  if (this%annotation%mesh /= '') then
328  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'mesh', &
329  this%annotation%mesh), this%nc_fname)
330 
331  end if
332  ! MODFLOW 6 model type
333  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'modflow_model', &
334  this%annotation%model), this%nc_fname)
335  ! generation datetime
336  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'history', &
337  this%annotation%history), this%nc_fname)
338  ! supported conventions
339  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'Conventions', &
340  this%annotation%conventions), &
341  this%nc_fname)
Here is the call graph for this function:

◆ add_pkg_data()

subroutine meshmodelmodule::add_pkg_data ( class(meshmodeltype), intent(inout)  this)

Definition at line 373 of file MeshNCModel.f90.

380  class(MeshModelType), intent(inout) :: this
381  character(LENCOMPONENTNAME) :: ptype, pname, pkgtype
382  character(len=LENMEMPATH) :: input_mempath
383  type(CharacterStringType), dimension(:), contiguous, &
384  pointer :: pkgtypes => null()
385  type(CharacterStringType), dimension(:), contiguous, &
386  pointer :: pkgnames => null()
387  type(CharacterStringType), dimension(:), contiguous, &
388  pointer :: mempaths => null()
389  type(InputParamDefinitionType), dimension(:), pointer :: param_dfns
390  character(len=LENMEMPATH) :: mempath
391  integer(I4B) :: n
392  integer(I4B), pointer :: export_arrays
393  logical(LGP) :: found
394 
395  input_mempath = create_mem_path(component=this%modelname, context=idm_context)
396 
397  ! set pointers to model path package info
398  call mem_setptr(pkgtypes, 'PKGTYPES', input_mempath)
399  call mem_setptr(pkgnames, 'PKGNAMES', input_mempath)
400  call mem_setptr(mempaths, 'MEMPATHS', input_mempath)
401 
402  allocate (export_arrays)
403 
404  do n = 1, size(mempaths)
405  ! initialize export_arrays
406  export_arrays = 0
407 
408  ! set package attributes
409  mempath = mempaths(n)
410  pname = pkgnames(n)
411  ptype = pkgtypes(n)
412 
413  ! export input arrays
414  if (mempath /= '') then
415  ! update export
416  call mem_set_value(export_arrays, 'EXPORT_NC', mempath, found)
417  if (export_arrays > 0) then
418  pkgtype = idm_subcomponent_type(this%modeltype, ptype)
419  param_dfns => param_definitions(this%modeltype, pkgtype)
420  call this%export_input_arrays(ptype, pname, mempath, param_dfns)
421  end if
422  end if
423  end do
424 
425  ! cleanup
426  deallocate (export_arrays)
type(inputparamdefinitiontype) function, dimension(:), pointer, public param_definitions(component, subcomponent)
logical function, public idm_multi_package(component, subcomponent)
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
character(len=lencomponentname) function, public idm_subcomponent_type(component, subcomponent)
component from package or model type
Here is the call graph for this function:

◆ create_mesh()

subroutine meshmodelmodule::create_mesh ( class(mesh2dmodeltype), intent(inout)  this)
private

Definition at line 516 of file MeshNCModel.f90.

517  class(Mesh2dModelType), intent(inout) :: this
518 
519  ! create mesh container variable
520  call nf_verify(nf90_def_var(this%ncid, this%mesh_name, nf90_int, &
521  this%var_ids%mesh), this%nc_fname)
522 
523  ! assign container variable attributes
524  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'cf_role', &
525  'mesh_topology'), this%nc_fname)
526  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'long_name', &
527  '2D mesh topology'), this%nc_fname)
528  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
529  'topology_dimension', 2), this%nc_fname)
530  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'face_dimension', &
531  'nmesh_face'), this%nc_fname)
532  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
533  'node_coordinates', 'mesh_node_x mesh_node_y'), &
534  this%nc_fname)
535  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
536  'face_coordinates', 'mesh_face_x mesh_face_y'), &
537  this%nc_fname)
538  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
539  'face_node_connectivity', 'mesh_face_nodes'), &
540  this%nc_fname)
541 
542  ! create mesh x node (mesh vertex) variable
543  call nf_verify(nf90_def_var(this%ncid, 'mesh_node_x', nf90_double, &
544  (/this%dim_ids%nmesh_node/), &
545  this%var_ids%mesh_node_x), this%nc_fname)
546 
547  ! assign mesh x node variable attributes
548  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
549  'units', this%lenunits), this%nc_fname)
550  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
551  'standard_name', 'projection_x_coordinate'), &
552  this%nc_fname)
553  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
554  'long_name', 'Easting'), this%nc_fname)
555 
556  if (this%wkt /= '') then
557  ! associate with projection
558  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
559  'grid_mapping', this%gridmap_name), &
560  this%nc_fname)
561  end if
562 
563  ! create mesh y node (mesh vertex) variable
564  call nf_verify(nf90_def_var(this%ncid, 'mesh_node_y', nf90_double, &
565  (/this%dim_ids%nmesh_node/), &
566  this%var_ids%mesh_node_y), this%nc_fname)
567 
568  ! assign mesh y variable attributes
569  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
570  'units', this%lenunits), this%nc_fname)
571  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
572  'standard_name', 'projection_y_coordinate'), &
573  this%nc_fname)
574  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
575  'long_name', 'Northing'), this%nc_fname)
576 
577  if (this%wkt /= '') then
578  ! associate with projection
579  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
580  'grid_mapping', this%gridmap_name), &
581  this%nc_fname)
582  end if
583 
584  ! create mesh x face (cell vertex) variable
585  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_x', nf90_double, &
586  (/this%dim_ids%nmesh_face/), &
587  this%var_ids%mesh_face_x), this%nc_fname)
588 
589  ! assign mesh x face variable attributes
590  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
591  'units', this%lenunits), this%nc_fname)
592  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
593  'standard_name', 'projection_x_coordinate'), &
594  this%nc_fname)
595  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
596  'long_name', 'Easting'), this%nc_fname)
597  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, 'bounds', &
598  'mesh_face_xbnds'), this%nc_fname)
599  if (this%wkt /= '') then
600  ! associate with projection
601  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
602  'grid_mapping', this%gridmap_name), &
603  this%nc_fname)
604  end if
605 
606  ! create mesh x cell bounds variable
607  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_xbnds', nf90_double, &
608  (/this%dim_ids%max_nmesh_face_nodes, &
609  this%dim_ids%nmesh_face/), &
610  this%var_ids%mesh_face_xbnds), &
611  this%nc_fname)
612 
613  ! create mesh y face (cell vertex) variable
614  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_y', nf90_double, &
615  (/this%dim_ids%nmesh_face/), &
616  this%var_ids%mesh_face_y), this%nc_fname)
617 
618  ! assign mesh y face variable attributes
619  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
620  'units', this%lenunits), this%nc_fname)
621  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
622  'standard_name', 'projection_y_coordinate'), &
623  this%nc_fname)
624  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
625  'long_name', 'Northing'), this%nc_fname)
626  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, 'bounds', &
627  'mesh_face_ybnds'), this%nc_fname)
628 
629  if (this%wkt /= '') then
630  ! associate with projection
631  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
632  'grid_mapping', this%gridmap_name), &
633  this%nc_fname)
634  end if
635 
636  ! create mesh y cell bounds variable
637  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_ybnds', nf90_double, &
638  (/this%dim_ids%max_nmesh_face_nodes, &
639  this%dim_ids%nmesh_face/), &
640  this%var_ids%mesh_face_ybnds), &
641  this%nc_fname)
642 
643  ! create mesh face nodes variable
644  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_nodes', nf90_int, &
645  (/this%dim_ids%max_nmesh_face_nodes, &
646  this%dim_ids%nmesh_face/), &
647  this%var_ids%mesh_face_nodes), &
648  this%nc_fname)
649 
650  ! assign variable attributes
651  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
652  'cf_role', 'face_node_connectivity'), &
653  this%nc_fname)
654  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
655  'long_name', &
656  'Vertices bounding cell (counterclockwise)'), &
657  this%nc_fname)
658  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
659  '_FillValue', (/nf90_fill_int/)), &
660  this%nc_fname)
661  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
662  'start_index', 1), this%nc_fname)
Here is the call graph for this function:

◆ create_timeseries()

subroutine meshmodelmodule::create_timeseries ( class(meshmodeltype), intent(inout)  this,
type(inputparamdefinitiontype), intent(in), pointer  idt,
integer(i4b), intent(in)  iparam,
integer(i4b), intent(in)  iaux,
integer(i4b), intent(in)  layer,
class(exportpackagetype), intent(in), pointer  export_pkg 
)

Definition at line 221 of file MeshNCModel.f90.

222  use constantsmodule, only: dnodata
224  class(MeshModelType), intent(inout) :: this
225  type(InputParamDefinitionType), pointer, intent(in) :: idt
226  integer(I4B), intent(in) :: iparam
227  integer(I4B), intent(in) :: iaux
228  integer(I4B), intent(in) :: layer
229  class(ExportPackageType), pointer, intent(in) :: export_pkg
230  character(len=LINELENGTH) :: varname, longname, nc_tag
231  integer(I4B) :: varid
232 
233  ! set variable input tag
234  nc_tag = this%input_attribute(export_pkg%mf6_input%subcomponent_name, &
235  idt)
236 
237  ! set names
238  varname = export_varname(export_pkg%mf6_input%subcomponent_name, &
239  idt%tagname, export_pkg%mf6_input%mempath, &
240  layer=layer, iaux=iaux)
241  longname = export_longname(idt%longname, &
242  export_pkg%mf6_input%subcomponent_name, &
243  idt%tagname, export_pkg%mf6_input%mempath, &
244  layer=layer, iaux=iaux)
245 
246  ! create the netcdf dependent layer variable
247  select case (idt%datatype)
248  case ('DOUBLE1D', 'DOUBLE2D')
249  call nf_verify(nf90_def_var(this%ncid, varname, nf90_double, &
250  (/this%dim_ids%nmesh_face, &
251  this%dim_ids%time/), &
252  varid), &
253  this%nc_fname)
254  call nf_verify(nf90_put_att(this%ncid, varid, &
255  '_FillValue', (/dnodata/)), &
256  this%nc_fname)
257  case ('INTEGER1D')
258  call nf_verify(nf90_def_var(this%ncid, varname, nf90_int, &
259  (/this%dim_ids%nmesh_face, &
260  this%dim_ids%time/), &
261  varid), &
262  this%nc_fname)
263  call nf_verify(nf90_put_att(this%ncid, varid, &
264  '_FillValue', (/nf90_fill_int/)), &
265  this%nc_fname)
266  end select
267 
268  ! apply chunking parameters
269  if (this%chunking_active) then
270  call nf_verify(nf90_def_var_chunking(this%ncid, &
271  varid, &
272  nf90_chunked, &
273  (/this%chunk_face, &
274  this%chunk_time/)), &
275  this%nc_fname)
276  end if
277 
278  ! deflate and shuffle
279  call ncvar_deflate(this%ncid, varid, this%deflate, &
280  this%shuffle, this%nc_fname)
281 
282  ! assign variable attributes
283  call nf_verify(nf90_put_att(this%ncid, varid, &
284  'units', this%lenunits), this%nc_fname)
285  call nf_verify(nf90_put_att(this%ncid, varid, &
286  'long_name', longname), this%nc_fname)
287  call nf_verify(nf90_put_att(this%ncid, varid, &
288  'mesh', this%mesh_name), this%nc_fname)
289  call nf_verify(nf90_put_att(this%ncid, varid, &
290  'location', 'face'), this%nc_fname)
291 
292  ! add grid mapping and mf6 attr
293  call ncvar_gridmap(this%ncid, varid, &
294  this%gridmap_name, this%nc_fname)
295  call ncvar_mf6attr(this%ncid, varid, layer, iaux, nc_tag, this%nc_fname)
296 
297  ! store variable id
298  if (idt%tagname == 'AUX') then
299  if (layer > 0) then
300  export_pkg%varids_aux(iaux, layer) = varid
301  else
302  export_pkg%varids_aux(iaux, 1) = varid
303  end if
304  else
305  if (layer > 0) then
306  export_pkg%varids_param(iparam, layer) = varid
307  else
308  export_pkg%varids_param(iparam, 1) = varid
309  end if
310  end if
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
This module contains the NCModelExportModule.
Definition: NCModel.f90:8
Here is the call graph for this function:

◆ define_dependent()

subroutine meshmodelmodule::define_dependent ( class(meshmodeltype), intent(inout)  this)

Definition at line 431 of file MeshNCModel.f90.

432  class(MeshModelType), intent(inout) :: this
433  character(len=LINELENGTH) :: varname, longname
434  integer(I4B) :: k
435 
436  ! create a dependent variable for each layer
437  do k = 1, this%nlay
438  ! initialize names
439  varname = ''
440  longname = ''
441 
442  ! set layer variable and longnames
443  write (varname, '(a,i0)') trim(this%xname)//'_l', k
444  write (longname, '(a,i0,a)') trim(this%annotation%longname)// &
445  ' (layer ', k, ')'
446 
447  ! create the netcdf dependent layer variable
448  call nf_verify(nf90_def_var(this%ncid, varname, nf90_double, &
449  (/this%dim_ids%nmesh_face, &
450  this%dim_ids%time/), &
451  this%var_ids%dependent(k)), &
452  this%nc_fname)
453 
454  ! apply chunking parameters
455  if (this%chunking_active) then
456  call nf_verify(nf90_def_var_chunking(this%ncid, &
457  this%var_ids%dependent(k), &
458  nf90_chunked, &
459  (/this%chunk_face, &
460  this%chunk_time/)), &
461  this%nc_fname)
462  end if
463 
464  ! deflate and shuffle
465  call ncvar_deflate(this%ncid, this%var_ids%dependent(k), this%deflate, &
466  this%shuffle, this%nc_fname)
467 
468  ! assign variable attributes
469  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
470  'units', this%lenunits), this%nc_fname)
471  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
472  'standard_name', this%annotation%stdname), &
473  this%nc_fname)
474  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
475  'long_name', longname), this%nc_fname)
476  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
477  '_FillValue', (/dhnoflo/)), &
478  this%nc_fname)
479  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
480  'mesh', this%mesh_name), this%nc_fname)
481  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
482  'location', 'face'), this%nc_fname)
483 
484  ! add grid mapping
485  call ncvar_gridmap(this%ncid, this%var_ids%dependent(k), &
486  this%gridmap_name, this%nc_fname)
487  end do
Here is the call graph for this function:

◆ define_gridmap()

subroutine meshmodelmodule::define_gridmap ( class(meshmodeltype), intent(inout)  this)
private

Definition at line 492 of file MeshNCModel.f90.

493  class(MeshModelType), intent(inout) :: this
494  integer(I4B) :: var_id
495 
496  ! was projection info provided
497  if (this%wkt /= '') then
498  ! create projection variable
499  call nf_verify(nf90_redef(this%ncid), this%nc_fname)
500  call nf_verify(nf90_def_var(this%ncid, this%gridmap_name, nf90_int, &
501  var_id), this%nc_fname)
502  ! cf-conventions prefers 'crs_wkt'
503  !call nf_verify(nf90_put_att(this%ncid, var_id, 'crs_wkt', this%wkt), &
504  ! this%nc_fname)
505  ! QGIS recognizes 'wkt'
506  call nf_verify(nf90_put_att(this%ncid, var_id, 'wkt', this%wkt), &
507  this%nc_fname)
508  call nf_verify(nf90_enddef(this%ncid), this%nc_fname)
509  call nf_verify(nf90_put_var(this%ncid, var_id, 1), &
510  this%nc_fname)
511  end if
Here is the call graph for this function:

◆ df_export()

subroutine meshmodelmodule::df_export ( class(meshmodeltype), intent(inout)  this)

Definition at line 164 of file MeshNCModel.f90.

166  class(MeshModelType), intent(inout) :: this
167  class(ExportPackageType), pointer :: export_pkg
168  integer(I4B) :: idx
169  do idx = 1, this%pkglist%Count()
170  export_pkg => this%get(idx)
171  call this%export_df(export_pkg)
172  end do

◆ export_df()

subroutine meshmodelmodule::export_df ( class(meshmodeltype), intent(inout)  this,
class(exportpackagetype), intent(in), pointer  export_pkg 
)

Definition at line 177 of file MeshNCModel.f90.

180  class(MeshModelType), intent(inout) :: this
181  class(ExportPackageType), pointer, intent(in) :: export_pkg
182  type(InputParamDefinitionType), pointer :: idt
183  integer(I4B) :: iparam, iaux, layer
184 
185  ! export defined period input
186  do iparam = 1, export_pkg%nparam
187  ! initialize
188  iaux = 0
189  layer = 0
190  ! set input definition
191  idt => &
192  get_param_definition_type(export_pkg%mf6_input%param_dfns, &
193  export_pkg%mf6_input%component_type, &
194  export_pkg%mf6_input%subcomponent_type, &
195  'PERIOD', export_pkg%param_names(iparam), '')
196 
197  select case (idt%shape)
198  case ('NCPL')
199  call this%create_timeseries(idt, iparam, iaux, layer, export_pkg)
200  case ('NODES')
201  do layer = 1, this%nlay
202  call this%create_timeseries(idt, iparam, iaux, layer, export_pkg)
203  end do
204  case ('NAUX NCPL')
205  do iaux = 1, export_pkg%naux
206  call this%create_timeseries(idt, iparam, iaux, layer, export_pkg)
207  end do
208  case ('NAUX NODES')
209  do iaux = 1, export_pkg%naux
210  do layer = 1, this%nlay
211  call this%create_timeseries(idt, iparam, iaux, layer, export_pkg)
212  end do
213  end do
214  case default
215  end select
216  end do
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.
Here is the call graph for this function:

◆ export_input_arrays()

subroutine meshmodelmodule::export_input_arrays ( class(meshmodeltype), intent(inout)  this,
character(len=*), intent(in)  pkgtype,
character(len=*), intent(in)  pkgname,
character(len=*), intent(in)  mempath,
type(inputparamdefinitiontype), dimension(:), intent(in), pointer  param_dfns 
)
private

Definition at line 346 of file MeshNCModel.f90.

347  use memorymanagermodule, only: get_isize
348  class(MeshModelType), intent(inout) :: this
349  character(len=*), intent(in) :: pkgtype
350  character(len=*), intent(in) :: pkgname
351  character(len=*), intent(in) :: mempath
352  type(InputParamDefinitionType), dimension(:), pointer, &
353  intent(in) :: param_dfns
354  type(InputParamDefinitionType), pointer :: idt
355  integer(I4B) :: iparam, isize
356  ! export griddata block parameters
357  do iparam = 1, size(param_dfns)
358  ! assign param definition pointer
359  idt => param_dfns(iparam)
360  ! for now only griddata is exported
361  if (idt%blockname == 'GRIDDATA') then
362  ! veriy variable is allocated
363  call get_isize(idt%mf6varname, mempath, isize)
364  if (isize > 0) then
365  call this%export_input_array(pkgtype, pkgname, mempath, idt)
366  end if
367  end if
368  end do
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
Here is the call graph for this function:

◆ mesh_destroy()

subroutine meshmodelmodule::mesh_destroy ( class(meshmodeltype), intent(inout)  this)

Definition at line 154 of file MeshNCModel.f90.

156  class(MeshModelType), intent(inout) :: this
157  call nf_verify(nf90_close(this%ncid), this%nc_fname)
158  deallocate (this%chunk_face)
159  nullify (this%chunk_face)
Here is the call graph for this function:

◆ mesh_init()

subroutine meshmodelmodule::mesh_init ( class(meshmodeltype), 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)  lenuni,
integer(i4b), intent(in)  iout 
)
private

Definition at line 102 of file MeshNCModel.f90.

105  class(MeshModelType), intent(inout) :: this
106  character(len=*), intent(in) :: modelname
107  character(len=*), intent(in) :: modeltype
108  character(len=*), intent(in) :: modelfname
109  character(len=*), intent(in) :: nc_fname
110  integer(I4B), intent(in) :: disenum
111  integer(I4B), intent(in) :: nctype
112  integer(I4B), intent(in) :: lenuni
113  integer(I4B), intent(in) :: iout
114  logical(LGP) :: found
115 
116  ! initialize base class
117  call this%NCModelExportType%init(modelname, modeltype, modelfname, nc_fname, &
118  disenum, nctype, iout)
119 
120  ! allocate and initialize
121  allocate (this%chunk_face)
122  this%chunk_face = -1
123 
124  ! update values from input context
125  if (this%ncf_mempath /= '') then
126  call mem_set_value(this%chunk_face, 'CHUNK_FACE', this%ncf_mempath, found)
127  end if
128 
129  if (this%chunk_time > 0 .and. this%chunk_face > 0) then
130  this%chunking_active = .true.
131  else if (this%chunk_time > 0 .or. this%chunk_face > 0) then
132  this%chunk_face = -1
133  this%chunk_time = -1
134  write (warnmsg, '(a)') 'Ignoring user provided NetCDF chunking parameter. &
135  &Define chunk_time and chunk_face input parameters to see an effect in &
136  &file "'//trim(nc_fname)//'".'
137  call store_warning(warnmsg)
138  end if
139 
140  if (lenuni == 1) then
141  this%lenunits = 'ft'
142  else
143  this%lenunits = 'm'
144  end if
145 
146  ! create the netcdf file
147  call nf_verify(nf90_create(this%nc_fname, &
148  ior(nf90_clobber, nf90_netcdf4), this%ncid), &
149  this%nc_fname)
Here is the call graph for this function:

◆ ncvar_chunk()

subroutine, public meshmodelmodule::ncvar_chunk ( integer(i4b), intent(in)  ncid,
integer(i4b), intent(in)  varid,
integer(i4b), intent(in)  chunk_face,
character(len=*), intent(in)  nc_fname 
)

Definition at line 667 of file MeshNCModel.f90.

668  integer(I4B), intent(in) :: ncid
669  integer(I4B), intent(in) :: varid
670  integer(I4B), intent(in) :: chunk_face
671  character(len=*), intent(in) :: nc_fname
672  if (chunk_face > 0) then
673  call nf_verify(nf90_def_var_chunking(ncid, varid, nf90_chunked, &
674  (/chunk_face/)), nc_fname)
675  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ncvar_deflate()

subroutine, public meshmodelmodule::ncvar_deflate ( integer(i4b), intent(in)  ncid,
integer(i4b), intent(in)  varid,
integer(i4b), intent(in)  deflate,
integer(i4b), intent(in)  shuffle,
character(len=*), intent(in)  nc_fname 
)

Definition at line 680 of file MeshNCModel.f90.

681  integer(I4B), intent(in) :: ncid
682  integer(I4B), intent(in) :: varid
683  integer(I4B), intent(in) :: deflate
684  integer(I4B), intent(in) :: shuffle
685  character(len=*), intent(in) :: nc_fname
686  if (deflate >= 0) then
687  call nf_verify(nf90_def_var_deflate(ncid, varid, shuffle=shuffle, &
688  deflate=1, deflate_level=deflate), &
689  nc_fname)
690  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ncvar_gridmap()

subroutine, public meshmodelmodule::ncvar_gridmap ( integer(i4b), intent(in)  ncid,
integer(i4b), intent(in)  varid,
character(len=*), intent(in)  gridmap_name,
character(len=*), intent(in)  nc_fname 
)

Definition at line 695 of file MeshNCModel.f90.

696  integer(I4B), intent(in) :: ncid
697  integer(I4B), intent(in) :: varid
698  character(len=*), intent(in) :: gridmap_name
699  character(len=*), intent(in) :: nc_fname
700  if (gridmap_name /= '') then
701  call nf_verify(nf90_put_att(ncid, varid, 'coordinates', &
702  'mesh_face_x mesh_face_y'), nc_fname)
703  call nf_verify(nf90_put_att(ncid, varid, 'grid_mapping', &
704  gridmap_name), nc_fname)
705  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ncvar_mf6attr()

subroutine, public meshmodelmodule::ncvar_mf6attr ( integer(i4b), intent(in)  ncid,
integer(i4b), intent(in)  varid,
integer(i4b), intent(in)  layer,
integer(i4b), intent(in)  iaux,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  nc_fname 
)

Definition at line 710 of file MeshNCModel.f90.

711  integer(I4B), intent(in) :: ncid
712  integer(I4B), intent(in) :: varid
713  integer(I4B), intent(in) :: layer
714  integer(I4B), intent(in) :: iaux
715  character(len=*), intent(in) :: nc_tag
716  character(len=*), intent(in) :: nc_fname
717  if (nc_tag /= '') then
718  call nf_verify(nf90_put_att(ncid, varid, 'modflow_input', &
719  nc_tag), nc_fname)
720  if (layer > 0) then
721  call nf_verify(nf90_put_att(ncid, varid, 'layer', &
722  layer), nc_fname)
723  end if
724  if (iaux > 0) then
725  call nf_verify(nf90_put_att(ncid, varid, 'modflow_iaux', &
726  iaux), nc_fname)
727  end if
728  end if
Here is the call graph for this function:
Here is the caller graph for this function: