MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
MeshNCModel.f90
Go to the documentation of this file.
1 !> @brief This module contains the MeshModelModule
2 !!
3 !! This module defines a base class for UGRID based
4 !! (mesh) model netcdf exports. It is dependent on
5 !! external netcdf libraries.
6 !<
8 
9  use kindmodule, only: dp, i4b, lgp
12  use simvariablesmodule, only: errmsg
18  use netcdfcommonmodule, only: nf_verify
19  use netcdf
20 
21  implicit none
22  private
24  public :: mesh2dmodeltype
25  public :: ncvar_chunk
26  public :: ncvar_deflate
27  public :: ncvar_gridmap
28  public :: ncvar_mf6attr
29  public :: export_varname
30 
31  !> @brief type for storing model export dimension ids
32  !<
34  integer(I4B) :: nmesh_node !< number of nodes in mesh
35  integer(I4B) :: nmesh_face !< number of faces in mesh
36  integer(I4B) :: max_nmesh_face_nodes !< max number of nodes in a single face
37  integer(I4B) :: nlay !< number of layers
38  integer(I4B) :: time !< number of steps
39  contains
40  end type meshncdimidtype
41 
42  !> @brief type for storing model export variable ids
43  !<
45  integer(I4B) :: mesh !< mesh container variable
46  integer(I4B) :: mesh_node_x !< mesh nodes x array
47  integer(I4B) :: mesh_node_y !< mesh nodes y array
48  integer(I4B) :: mesh_face_x !< mesh faces x location array
49  integer(I4B) :: mesh_face_y !< mesh faces y location array
50  integer(I4B) :: mesh_face_xbnds !< mesh faces 2D x bounds array
51  integer(I4B) :: mesh_face_ybnds !< mesh faces 2D y bounds array
52  integer(I4B) :: mesh_face_nodes !< mesh faces 2D nodes array
53  integer(I4B) :: time !< time coordinate variable
54  integer(I4B), dimension(:), allocatable :: dependent !< layered dependent variables array
55  contains
56  end type meshncvaridtype
57 
58  !> @brief base ugrid netcdf export type
59  !<
60  type, abstract, extends(ncbasemodelexporttype) :: meshmodeltype
61  type(meshncdimidtype) :: dim_ids !< dimension ids
62  type(meshncvaridtype) :: var_ids !< variable ids
63  integer(I4B) :: nlay !< number of layers
64  integer(I4B), pointer :: chunk_face !< chunking parameter for face dimension
65  contains
66  procedure :: mesh_init
67  procedure :: mesh_destroy
68  procedure :: add_global_att
69  procedure(nc_array_export_if), deferred :: export_input_array
70  procedure :: export_input_arrays
71  procedure :: add_pkg_data
72  procedure :: define_dependent
73  procedure :: define_gridmap
74  end type meshmodeltype
75 
76  !> @brief abstract interfaces for derived ugrid netcd export types
77  !<
78  abstract interface
79  subroutine nc_array_export_if(this, pkgtype, pkgname, mempath, idt)
81  class(meshmodeltype), intent(inout) :: this
82  character(len=*), intent(in) :: pkgtype
83  character(len=*), intent(in) :: pkgname
84  character(len=*), intent(in) :: mempath
85  type(inputparamdefinitiontype), pointer, intent(in) :: idt
86  end subroutine
87  end interface
88 
89  type, abstract, extends(meshmodeltype) :: mesh2dmodeltype
90  contains
91  procedure :: create_mesh
92  end type mesh2dmodeltype
93 
94 contains
95 
96  !> @brief initialize
97  !<
98  subroutine mesh_init(this, modelname, modeltype, modelfname, nc_fname, &
99  disenum, nctype, iout)
101  class(meshmodeltype), intent(inout) :: this
102  character(len=*), intent(in) :: modelname
103  character(len=*), intent(in) :: modeltype
104  character(len=*), intent(in) :: modelfname
105  character(len=*), intent(in) :: nc_fname
106  integer(I4B), intent(in) :: disenum
107  integer(I4B), intent(in) :: nctype
108  integer(I4B), intent(in) :: iout
109  logical(LGP) :: found
110 
111  ! initialize base class
112  call this%NCModelExportType%init(modelname, modeltype, modelfname, nc_fname, &
113  disenum, nctype, iout)
114 
115  ! allocate and initialize
116  allocate (this%chunk_face)
117  this%chunk_face = -1
118 
119  ! update values from input context
120  if (this%ncf_mempath /= '') then
121  call mem_set_value(this%chunk_face, 'CHUNK_FACE', this%ncf_mempath, found)
122  end if
123 
124  if (this%chunk_time > 0 .and. this%chunk_face > 0) then
125  this%chunking_active = .true.
126  end if
127 
128  ! create the netcdf file
129  call nf_verify(nf90_create(this%nc_fname, &
130  iand(nf90_clobber, nf90_netcdf4), this%ncid), &
131  this%nc_fname)
132  end subroutine mesh_init
133 
134  !> @brief initialize
135  !<
136  subroutine mesh_destroy(this)
138  class(meshmodeltype), intent(inout) :: this
139  call nf_verify(nf90_close(this%ncid), this%nc_fname)
140  deallocate (this%chunk_face)
141  nullify (this%chunk_face)
142  end subroutine mesh_destroy
143 
144  !> @brief create file (group) attributes
145  !<
146  subroutine add_global_att(this)
147  class(meshmodeltype), intent(inout) :: this
148  ! file scoped title
149  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'title', &
150  this%annotation%title), this%nc_fname)
151  ! source (MODFLOW 6)
152  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'source', &
153  this%annotation%source), this%nc_fname)
154  ! export type (MODFLOW 6)
155  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'modflow6_grid', &
156  this%annotation%grid), this%nc_fname)
157  ! MODFLOW 6 model type
158  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'modflow6_model', &
159  this%annotation%model), this%nc_fname)
160  ! generation datetime
161  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'history', &
162  this%annotation%history), this%nc_fname)
163  ! supported conventions
164  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'Conventions', &
165  this%annotation%conventions), &
166  this%nc_fname)
167  end subroutine add_global_att
168 
169  !> @brief write package gridded input data
170  !<
171  subroutine export_input_arrays(this, pkgtype, pkgname, mempath, param_dfns)
172  use memorymanagermodule, only: get_isize
173  class(meshmodeltype), intent(inout) :: this
174  character(len=*), intent(in) :: pkgtype
175  character(len=*), intent(in) :: pkgname
176  character(len=*), intent(in) :: mempath
177  type(inputparamdefinitiontype), dimension(:), pointer, &
178  intent(in) :: param_dfns
179  type(inputparamdefinitiontype), pointer :: idt
180  integer(I4B) :: iparam, isize
181  ! export griddata block parameters
182  do iparam = 1, size(param_dfns)
183  ! assign param definition pointer
184  idt => param_dfns(iparam)
185  ! for now only griddata is exported
186  if (idt%blockname == 'GRIDDATA') then
187  ! veriy variable is allocated
188  call get_isize(idt%mf6varname, mempath, isize)
189  if (isize > 0) then
190  call this%export_input_array(pkgtype, pkgname, mempath, idt)
191  end if
192  end if
193  end do
194  end subroutine export_input_arrays
195 
196  !> @brief determine packages to write gridded input
197  !<
198  subroutine add_pkg_data(this)
205  class(meshmodeltype), intent(inout) :: this
206  character(LENCOMPONENTNAME) :: ptype, pname, pkgtype
207  character(len=LENMEMPATH) :: input_mempath
208  type(characterstringtype), dimension(:), contiguous, &
209  pointer :: pkgtypes => null()
210  type(characterstringtype), dimension(:), contiguous, &
211  pointer :: pkgnames => null()
212  type(characterstringtype), dimension(:), contiguous, &
213  pointer :: mempaths => null()
214  type(inputparamdefinitiontype), dimension(:), pointer :: param_dfns
215  character(len=LENMEMPATH) :: mempath
216  integer(I4B) :: n
217  integer(I4B), pointer :: export_arrays
218  logical(LGP) :: found
219 
220  input_mempath = create_mem_path(component=this%modelname, context=idm_context)
221 
222  ! set pointers to model path package info
223  call mem_setptr(pkgtypes, 'PKGTYPES', input_mempath)
224  call mem_setptr(pkgnames, 'PKGNAMES', input_mempath)
225  call mem_setptr(mempaths, 'MEMPATHS', input_mempath)
226 
227  allocate (export_arrays)
228 
229  do n = 1, size(mempaths)
230  ! initialize export_arrays
231  export_arrays = 0
232 
233  ! set package attributes
234  mempath = mempaths(n)
235  pname = pkgnames(n)
236  ptype = pkgtypes(n)
237 
238  ! export input arrays
239  if (mempath /= '') then
240  ! update export
241  call mem_set_value(export_arrays, 'EXPORT_NC', mempath, found)
242  if (export_arrays > 0) then
243  pkgtype = idm_subcomponent_type(this%modeltype, ptype)
244  param_dfns => param_definitions(this%modeltype, pkgtype)
245  call this%export_input_arrays(ptype, pname, mempath, param_dfns)
246  end if
247  end if
248  end do
249 
250  ! cleanup
251  deallocate (export_arrays)
252  end subroutine add_pkg_data
253 
254  !> @brief create the model layer dependent variables
255  !<
256  subroutine define_dependent(this)
257  class(meshmodeltype), intent(inout) :: this
258  character(len=LINELENGTH) :: varname, longname
259  integer(I4B) :: k
260 
261  ! create a dependent variable for each layer
262  do k = 1, this%nlay
263  ! initialize names
264  varname = ''
265  longname = ''
266 
267  ! set layer variable and longnames
268  write (varname, '(a,i0)') trim(this%xname)//'_l', k
269  write (longname, '(a,i0,a)') trim(this%annotation%longname)// &
270  ' (layer ', k, ')'
271 
272  ! create the netcdf dependent layer variable
273  call nf_verify(nf90_def_var(this%ncid, varname, nf90_double, &
274  (/this%dim_ids%nmesh_face, &
275  this%dim_ids%time/), &
276  this%var_ids%dependent(k)), &
277  this%nc_fname)
278 
279  ! apply chunking parameters
280  if (this%chunking_active) then
281  call nf_verify(nf90_def_var_chunking(this%ncid, &
282  this%var_ids%dependent(k), &
283  nf90_chunked, &
284  (/this%chunk_face, &
285  this%chunk_time/)), &
286  this%nc_fname)
287  end if
288 
289  ! deflate and shuffle
290  call ncvar_deflate(this%ncid, this%var_ids%dependent(k), this%deflate, &
291  this%shuffle, this%nc_fname)
292 
293  ! assign variable attributes
294  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
295  'units', 'm'), this%nc_fname)
296  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
297  'standard_name', this%annotation%stdname), &
298  this%nc_fname)
299  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
300  'long_name', longname), this%nc_fname)
301  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
302  '_FillValue', (/dhnoflo/)), &
303  this%nc_fname)
304  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
305  'mesh', this%mesh_name), this%nc_fname)
306  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
307  'location', 'face'), this%nc_fname)
308 
309  ! add grid mapping
310  call ncvar_gridmap(this%ncid, this%var_ids%dependent(k), &
311  this%gridmap_name, this%nc_fname)
312  end do
313  end subroutine define_dependent
314 
315  !> @brief create the file grid mapping container variable
316  !<
317  subroutine define_gridmap(this)
318  class(meshmodeltype), intent(inout) :: this
319  integer(I4B) :: var_id
320 
321  ! was projection info provided
322  if (this%ogc_wkt /= '') then
323  ! create projection variable
324  call nf_verify(nf90_redef(this%ncid), this%nc_fname)
325  call nf_verify(nf90_def_var(this%ncid, this%gridmap_name, nf90_int, &
326  var_id), this%nc_fname)
327  ! cf-conventions prefers 'crs_wkt'
328  !call nf_verify(nf90_put_att(this%ncid, var_id, 'crs_wkt', this%ogc_wkt), &
329  ! this%nc_fname)
330  ! QGIS recognizes 'wkt'
331  call nf_verify(nf90_put_att(this%ncid, var_id, 'wkt', this%ogc_wkt), &
332  this%nc_fname)
333  call nf_verify(nf90_enddef(this%ncid), this%nc_fname)
334  call nf_verify(nf90_put_var(this%ncid, var_id, 1), &
335  this%nc_fname)
336  end if
337  end subroutine define_gridmap
338 
339  !> @brief create the file mesh container variable
340  !<
341  subroutine create_mesh(this)
342  class(mesh2dmodeltype), intent(inout) :: this
343 
344  ! create mesh container variable
345  call nf_verify(nf90_def_var(this%ncid, this%mesh_name, nf90_int, &
346  this%var_ids%mesh), this%nc_fname)
347 
348  ! assign container variable attributes
349  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'cf_role', &
350  'mesh_topology'), this%nc_fname)
351  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'long_name', &
352  '2D mesh topology'), this%nc_fname)
353  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
354  'topology_dimension', 2), this%nc_fname)
355  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'face_dimension', &
356  'nmesh_face'), this%nc_fname)
357  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
358  'node_coordinates', 'mesh_node_x mesh_node_y'), &
359  this%nc_fname)
360  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
361  'face_coordinates', 'mesh_face_x mesh_face_y'), &
362  this%nc_fname)
363  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
364  'face_node_connectivity', 'mesh_face_nodes'), &
365  this%nc_fname)
366 
367  ! create mesh x node (mesh vertex) variable
368  call nf_verify(nf90_def_var(this%ncid, 'mesh_node_x', nf90_double, &
369  (/this%dim_ids%nmesh_node/), &
370  this%var_ids%mesh_node_x), this%nc_fname)
371 
372  ! assign mesh x node variable attributes
373  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
374  'units', 'm'), this%nc_fname)
375  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
376  'standard_name', 'projection_x_coordinate'), &
377  this%nc_fname)
378  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
379  'long_name', 'Easting'), this%nc_fname)
380 
381  if (this%ogc_wkt /= '') then
382  ! associate with projection
383  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
384  'grid_mapping', this%gridmap_name), &
385  this%nc_fname)
386  end if
387 
388  ! create mesh y node (mesh vertex) variable
389  call nf_verify(nf90_def_var(this%ncid, 'mesh_node_y', nf90_double, &
390  (/this%dim_ids%nmesh_node/), &
391  this%var_ids%mesh_node_y), this%nc_fname)
392 
393  ! assign mesh y variable attributes
394  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
395  'units', 'm'), this%nc_fname)
396  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
397  'standard_name', 'projection_y_coordinate'), &
398  this%nc_fname)
399  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
400  'long_name', 'Northing'), this%nc_fname)
401 
402  if (this%ogc_wkt /= '') then
403  ! associate with projection
404  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
405  'grid_mapping', this%gridmap_name), &
406  this%nc_fname)
407  end if
408 
409  ! create mesh x face (cell vertex) variable
410  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_x', nf90_double, &
411  (/this%dim_ids%nmesh_face/), &
412  this%var_ids%mesh_face_x), this%nc_fname)
413 
414  ! assign mesh x face variable attributes
415  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
416  'units', 'm'), this%nc_fname)
417  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
418  'standard_name', 'projection_x_coordinate'), &
419  this%nc_fname)
420  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
421  'long_name', 'Easting'), this%nc_fname)
422  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, 'bounds', &
423  'mesh_face_xbnds'), this%nc_fname)
424  if (this%ogc_wkt /= '') then
425  ! associate with projection
426  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
427  'grid_mapping', this%gridmap_name), &
428  this%nc_fname)
429  end if
430 
431  ! create mesh x cell bounds variable
432  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_xbnds', nf90_double, &
433  (/this%dim_ids%max_nmesh_face_nodes, &
434  this%dim_ids%nmesh_face/), &
435  this%var_ids%mesh_face_xbnds), &
436  this%nc_fname)
437 
438  ! create mesh y face (cell vertex) variable
439  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_y', nf90_double, &
440  (/this%dim_ids%nmesh_face/), &
441  this%var_ids%mesh_face_y), this%nc_fname)
442 
443  ! assign mesh y face variable attributes
444  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
445  'units', 'm'), this%nc_fname)
446  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
447  'standard_name', 'projection_y_coordinate'), &
448  this%nc_fname)
449  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
450  'long_name', 'Northing'), this%nc_fname)
451  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, 'bounds', &
452  'mesh_face_ybnds'), this%nc_fname)
453 
454  if (this%ogc_wkt /= '') then
455  ! associate with projection
456  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
457  'grid_mapping', this%gridmap_name), &
458  this%nc_fname)
459  end if
460 
461  ! create mesh y cell bounds variable
462  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_ybnds', nf90_double, &
463  (/this%dim_ids%max_nmesh_face_nodes, &
464  this%dim_ids%nmesh_face/), &
465  this%var_ids%mesh_face_ybnds), &
466  this%nc_fname)
467 
468  ! create mesh face nodes variable
469  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_nodes', nf90_int, &
470  (/this%dim_ids%max_nmesh_face_nodes, &
471  this%dim_ids%nmesh_face/), &
472  this%var_ids%mesh_face_nodes), &
473  this%nc_fname)
474 
475  ! assign variable attributes
476  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
477  'cf_role', 'face_node_connectivity'), &
478  this%nc_fname)
479  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
480  'long_name', &
481  'Vertices bounding cell (counterclockwise)'), &
482  this%nc_fname)
483  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
484  '_FillValue', (/nf90_fill_int/)), &
485  this%nc_fname)
486  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
487  'start_index', 1), this%nc_fname)
488  end subroutine create_mesh
489 
490  !> @brief define variable chunking
491  !<
492  subroutine ncvar_chunk(ncid, varid, chunk_face, nc_fname)
493  integer(I4B), intent(in) :: ncid
494  integer(I4B), intent(in) :: varid
495  integer(I4B), intent(in) :: chunk_face
496  character(len=*), intent(in) :: nc_fname
497  if (chunk_face > 0) then
498  call nf_verify(nf90_def_var_chunking(ncid, varid, nf90_chunked, &
499  (/chunk_face/)), nc_fname)
500  end if
501  end subroutine ncvar_chunk
502 
503  !> @brief define variable compression
504  !<
505  subroutine ncvar_deflate(ncid, varid, deflate, shuffle, nc_fname)
506  integer(I4B), intent(in) :: ncid
507  integer(I4B), intent(in) :: varid
508  integer(I4B), intent(in) :: deflate
509  integer(I4B), intent(in) :: shuffle
510  character(len=*), intent(in) :: nc_fname
511  if (deflate >= 0) then
512  call nf_verify(nf90_def_var_deflate(ncid, varid, shuffle=shuffle, &
513  deflate=1, deflate_level=deflate), &
514  nc_fname)
515  end if
516  end subroutine ncvar_deflate
517 
518  !> @brief put variable gridmap attributes
519  !<
520  subroutine ncvar_gridmap(ncid, varid, gridmap_name, nc_fname)
521  integer(I4B), intent(in) :: ncid
522  integer(I4B), intent(in) :: varid
523  character(len=*), intent(in) :: gridmap_name
524  character(len=*), intent(in) :: nc_fname
525  if (gridmap_name /= '') then
526  call nf_verify(nf90_put_att(ncid, varid, 'coordinates', &
527  'mesh_face_x mesh_face_y'), nc_fname)
528  call nf_verify(nf90_put_att(ncid, varid, 'grid_mapping', &
529  gridmap_name), nc_fname)
530  end if
531  end subroutine ncvar_gridmap
532 
533  !> @brief put variable internal attributes
534  !<
535  subroutine ncvar_mf6attr(ncid, varid, layer, iper, iaux, nc_tag, nc_fname)
536  integer(I4B), intent(in) :: ncid
537  integer(I4B), intent(in) :: varid
538  integer(I4B), intent(in) :: layer
539  integer(I4B), intent(in) :: iper
540  integer(I4B), intent(in) :: iaux
541  character(len=*), intent(in) :: nc_tag
542  character(len=*), intent(in) :: nc_fname
543  if (nc_tag /= '') then
544  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_input', &
545  nc_tag), nc_fname)
546  if (layer > 0) then
547  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_layer', &
548  layer), nc_fname)
549  end if
550  if (iper > 0) then
551  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_iper', &
552  iper), nc_fname)
553  end if
554  if (iaux > 0) then
555  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_iaux', &
556  iaux), nc_fname)
557  end if
558  end if
559  end subroutine ncvar_mf6attr
560 
561  !> @brief build netcdf variable name
562  !<
563  function export_varname(varname, layer, iper, iaux) result(vname)
564  use inputoutputmodule, only: lowcase
565  character(len=*), intent(in) :: varname
566  integer(I4B), optional, intent(in) :: layer
567  integer(I4B), optional, intent(in) :: iper
568  integer(I4B), optional, intent(in) :: iaux
569  character(len=LINELENGTH) :: vname
570  vname = ''
571  if (varname /= '') then
572  vname = varname
573  call lowcase(vname)
574  if (present(layer)) then
575  if (layer > 0) then
576  write (vname, '(a,i0)') trim(vname)//'_l', layer
577  end if
578  end if
579  if (present(iper)) then
580  if (iper > 0) then
581  write (vname, '(a,i0)') trim(vname)//'_p', iper
582  end if
583  end if
584  if (present(iaux)) then
585  if (iaux > 0) then
586  write (vname, '(a,i0)') trim(vname)//'a', iaux
587  end if
588  end if
589  end if
590  end function export_varname
591 
592 end module meshmodelmodule
abstract interfaces for derived ugrid netcd export types
Definition: MeshNCModel.f90:79
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
integer(i4b), parameter lencomponentname
maximum length of a component name
Definition: Constants.f90:18
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
real(dp), parameter dhnoflo
real no flow constant
Definition: Constants.f90:93
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
type(inputparamdefinitiontype) function, dimension(:), pointer, public param_definitions(component, subcomponent)
logical function, public idm_multi_package(component, subcomponent)
This module contains the InputDefinitionModule.
subroutine, public lowcase(word)
Convert to lower case.
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
subroutine, public get_isize(name, mem_path, isize)
@ brief Get the number of elements for this variable
This module contains the MeshModelModule.
Definition: MeshNCModel.f90:7
subroutine define_gridmap(this)
create the file grid mapping container variable
character(len=linelength) function, public export_varname(varname, layer, iper, iaux)
build netcdf variable name
subroutine, public ncvar_gridmap(ncid, varid, gridmap_name, nc_fname)
put variable gridmap attributes
subroutine mesh_init(this, modelname, modeltype, modelfname, nc_fname, disenum, nctype, iout)
initialize
subroutine, public ncvar_chunk(ncid, varid, chunk_face, nc_fname)
define variable chunking
subroutine mesh_destroy(this)
initialize
subroutine add_global_att(this)
create file (group) attributes
subroutine add_pkg_data(this)
determine packages to write gridded input
subroutine export_input_arrays(this, pkgtype, pkgname, mempath, param_dfns)
write package gridded input data
subroutine define_dependent(this)
create the model layer dependent variables
subroutine, public ncvar_deflate(ncid, varid, deflate, shuffle, nc_fname)
define variable compression
subroutine, public ncvar_mf6attr(ncid, varid, layer, iper, iaux, nc_tag, nc_fname)
put variable internal attributes
subroutine create_mesh(this)
create the file mesh container variable
This module contains the NCModelExportModule.
Definition: NCModel.f90:8
This module contains the NetCDFCommonModule.
Definition: NetCDFCommon.f90:6
subroutine, public nf_verify(res, nc_fname)
error check a netcdf-fortran interface call
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
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
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
base ugrid netcdf export type
Definition: MeshNCModel.f90:60
type for storing model export dimension ids
Definition: MeshNCModel.f90:33
type for storing model export variable ids
Definition: MeshNCModel.f90:44
abstract type for model netcdf export type
Definition: NCModel.f90:101