MODFLOW 6  version 6.7.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
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  else if (this%chunk_time > 0 .or. this%chunk_face > 0) then
127  this%chunk_face = -1
128  this%chunk_time = -1
129  write (warnmsg, '(a)') 'Ignoring user provided NetCDF chunking parameter. &
130  &Define chunk_time and chunk_face input parameters to see an effect.'
131  call store_warning(warnmsg)
132  end if
133 
134  ! create the netcdf file
135  call nf_verify(nf90_create(this%nc_fname, &
136  ior(nf90_clobber, nf90_netcdf4), this%ncid), &
137  this%nc_fname)
138  end subroutine mesh_init
139 
140  !> @brief initialize
141  !<
142  subroutine mesh_destroy(this)
144  class(meshmodeltype), intent(inout) :: this
145  call nf_verify(nf90_close(this%ncid), this%nc_fname)
146  deallocate (this%chunk_face)
147  nullify (this%chunk_face)
148  end subroutine mesh_destroy
149 
150  !> @brief create file (group) attributes
151  !<
152  subroutine add_global_att(this)
153  class(meshmodeltype), intent(inout) :: this
154  ! file scoped title
155  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'title', &
156  this%annotation%title), this%nc_fname)
157  ! source (MODFLOW 6)
158  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'source', &
159  this%annotation%source), this%nc_fname)
160  ! export type (MODFLOW 6)
161  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'modflow6_grid', &
162  this%annotation%grid), this%nc_fname)
163  ! MODFLOW 6 model type
164  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'modflow6_model', &
165  this%annotation%model), this%nc_fname)
166  ! generation datetime
167  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'history', &
168  this%annotation%history), this%nc_fname)
169  ! supported conventions
170  call nf_verify(nf90_put_att(this%ncid, nf90_global, 'Conventions', &
171  this%annotation%conventions), &
172  this%nc_fname)
173  end subroutine add_global_att
174 
175  !> @brief write package gridded input data
176  !<
177  subroutine export_input_arrays(this, pkgtype, pkgname, mempath, param_dfns)
178  use memorymanagermodule, only: get_isize
179  class(meshmodeltype), intent(inout) :: this
180  character(len=*), intent(in) :: pkgtype
181  character(len=*), intent(in) :: pkgname
182  character(len=*), intent(in) :: mempath
183  type(inputparamdefinitiontype), dimension(:), pointer, &
184  intent(in) :: param_dfns
185  type(inputparamdefinitiontype), pointer :: idt
186  integer(I4B) :: iparam, isize
187  ! export griddata block parameters
188  do iparam = 1, size(param_dfns)
189  ! assign param definition pointer
190  idt => param_dfns(iparam)
191  ! for now only griddata is exported
192  if (idt%blockname == 'GRIDDATA') then
193  ! veriy variable is allocated
194  call get_isize(idt%mf6varname, mempath, isize)
195  if (isize > 0) then
196  call this%export_input_array(pkgtype, pkgname, mempath, idt)
197  end if
198  end if
199  end do
200  end subroutine export_input_arrays
201 
202  !> @brief determine packages to write gridded input
203  !<
204  subroutine add_pkg_data(this)
211  class(meshmodeltype), intent(inout) :: this
212  character(LENCOMPONENTNAME) :: ptype, pname, pkgtype
213  character(len=LENMEMPATH) :: input_mempath
214  type(characterstringtype), dimension(:), contiguous, &
215  pointer :: pkgtypes => null()
216  type(characterstringtype), dimension(:), contiguous, &
217  pointer :: pkgnames => null()
218  type(characterstringtype), dimension(:), contiguous, &
219  pointer :: mempaths => null()
220  type(inputparamdefinitiontype), dimension(:), pointer :: param_dfns
221  character(len=LENMEMPATH) :: mempath
222  integer(I4B) :: n
223  integer(I4B), pointer :: export_arrays
224  logical(LGP) :: found
225 
226  input_mempath = create_mem_path(component=this%modelname, context=idm_context)
227 
228  ! set pointers to model path package info
229  call mem_setptr(pkgtypes, 'PKGTYPES', input_mempath)
230  call mem_setptr(pkgnames, 'PKGNAMES', input_mempath)
231  call mem_setptr(mempaths, 'MEMPATHS', input_mempath)
232 
233  allocate (export_arrays)
234 
235  do n = 1, size(mempaths)
236  ! initialize export_arrays
237  export_arrays = 0
238 
239  ! set package attributes
240  mempath = mempaths(n)
241  pname = pkgnames(n)
242  ptype = pkgtypes(n)
243 
244  ! export input arrays
245  if (mempath /= '') then
246  ! update export
247  call mem_set_value(export_arrays, 'EXPORT_NC', mempath, found)
248  if (export_arrays > 0) then
249  pkgtype = idm_subcomponent_type(this%modeltype, ptype)
250  param_dfns => param_definitions(this%modeltype, pkgtype)
251  call this%export_input_arrays(ptype, pname, mempath, param_dfns)
252  end if
253  end if
254  end do
255 
256  ! cleanup
257  deallocate (export_arrays)
258  end subroutine add_pkg_data
259 
260  !> @brief create the model layer dependent variables
261  !<
262  subroutine define_dependent(this)
263  class(meshmodeltype), intent(inout) :: this
264  character(len=LINELENGTH) :: varname, longname
265  integer(I4B) :: k
266 
267  ! create a dependent variable for each layer
268  do k = 1, this%nlay
269  ! initialize names
270  varname = ''
271  longname = ''
272 
273  ! set layer variable and longnames
274  write (varname, '(a,i0)') trim(this%xname)//'_l', k
275  write (longname, '(a,i0,a)') trim(this%annotation%longname)// &
276  ' (layer ', k, ')'
277 
278  ! create the netcdf dependent layer variable
279  call nf_verify(nf90_def_var(this%ncid, varname, nf90_double, &
280  (/this%dim_ids%nmesh_face, &
281  this%dim_ids%time/), &
282  this%var_ids%dependent(k)), &
283  this%nc_fname)
284 
285  ! apply chunking parameters
286  if (this%chunking_active) then
287  call nf_verify(nf90_def_var_chunking(this%ncid, &
288  this%var_ids%dependent(k), &
289  nf90_chunked, &
290  (/this%chunk_face, &
291  this%chunk_time/)), &
292  this%nc_fname)
293  end if
294 
295  ! deflate and shuffle
296  call ncvar_deflate(this%ncid, this%var_ids%dependent(k), this%deflate, &
297  this%shuffle, this%nc_fname)
298 
299  ! assign variable attributes
300  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
301  'units', 'm'), this%nc_fname)
302  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
303  'standard_name', this%annotation%stdname), &
304  this%nc_fname)
305  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
306  'long_name', longname), this%nc_fname)
307  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
308  '_FillValue', (/dhnoflo/)), &
309  this%nc_fname)
310  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
311  'mesh', this%mesh_name), this%nc_fname)
312  call nf_verify(nf90_put_att(this%ncid, this%var_ids%dependent(k), &
313  'location', 'face'), this%nc_fname)
314 
315  ! add grid mapping
316  call ncvar_gridmap(this%ncid, this%var_ids%dependent(k), &
317  this%gridmap_name, this%nc_fname)
318  end do
319  end subroutine define_dependent
320 
321  !> @brief create the file grid mapping container variable
322  !<
323  subroutine define_gridmap(this)
324  class(meshmodeltype), intent(inout) :: this
325  integer(I4B) :: var_id
326 
327  ! was projection info provided
328  if (this%wkt /= '') then
329  ! create projection variable
330  call nf_verify(nf90_redef(this%ncid), this%nc_fname)
331  call nf_verify(nf90_def_var(this%ncid, this%gridmap_name, nf90_int, &
332  var_id), this%nc_fname)
333  ! cf-conventions prefers 'crs_wkt'
334  !call nf_verify(nf90_put_att(this%ncid, var_id, 'crs_wkt', this%wkt), &
335  ! this%nc_fname)
336  ! QGIS recognizes 'wkt'
337  call nf_verify(nf90_put_att(this%ncid, var_id, 'wkt', this%wkt), &
338  this%nc_fname)
339  call nf_verify(nf90_enddef(this%ncid), this%nc_fname)
340  call nf_verify(nf90_put_var(this%ncid, var_id, 1), &
341  this%nc_fname)
342  end if
343  end subroutine define_gridmap
344 
345  !> @brief create the file mesh container variable
346  !<
347  subroutine create_mesh(this)
348  class(mesh2dmodeltype), intent(inout) :: this
349 
350  ! create mesh container variable
351  call nf_verify(nf90_def_var(this%ncid, this%mesh_name, nf90_int, &
352  this%var_ids%mesh), this%nc_fname)
353 
354  ! assign container variable attributes
355  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'cf_role', &
356  'mesh_topology'), this%nc_fname)
357  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'long_name', &
358  '2D mesh topology'), this%nc_fname)
359  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
360  'topology_dimension', 2), this%nc_fname)
361  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, 'face_dimension', &
362  'nmesh_face'), this%nc_fname)
363  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
364  'node_coordinates', 'mesh_node_x mesh_node_y'), &
365  this%nc_fname)
366  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
367  'face_coordinates', 'mesh_face_x mesh_face_y'), &
368  this%nc_fname)
369  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh, &
370  'face_node_connectivity', 'mesh_face_nodes'), &
371  this%nc_fname)
372 
373  ! create mesh x node (mesh vertex) variable
374  call nf_verify(nf90_def_var(this%ncid, 'mesh_node_x', nf90_double, &
375  (/this%dim_ids%nmesh_node/), &
376  this%var_ids%mesh_node_x), this%nc_fname)
377 
378  ! assign mesh x node variable attributes
379  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
380  'units', 'm'), this%nc_fname)
381  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
382  'standard_name', 'projection_x_coordinate'), &
383  this%nc_fname)
384  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
385  'long_name', 'Easting'), this%nc_fname)
386 
387  if (this%wkt /= '') then
388  ! associate with projection
389  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_x, &
390  'grid_mapping', this%gridmap_name), &
391  this%nc_fname)
392  end if
393 
394  ! create mesh y node (mesh vertex) variable
395  call nf_verify(nf90_def_var(this%ncid, 'mesh_node_y', nf90_double, &
396  (/this%dim_ids%nmesh_node/), &
397  this%var_ids%mesh_node_y), this%nc_fname)
398 
399  ! assign mesh y variable attributes
400  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
401  'units', 'm'), this%nc_fname)
402  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
403  'standard_name', 'projection_y_coordinate'), &
404  this%nc_fname)
405  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
406  'long_name', 'Northing'), this%nc_fname)
407 
408  if (this%wkt /= '') then
409  ! associate with projection
410  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_node_y, &
411  'grid_mapping', this%gridmap_name), &
412  this%nc_fname)
413  end if
414 
415  ! create mesh x face (cell vertex) variable
416  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_x', nf90_double, &
417  (/this%dim_ids%nmesh_face/), &
418  this%var_ids%mesh_face_x), this%nc_fname)
419 
420  ! assign mesh x face variable attributes
421  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
422  'units', 'm'), this%nc_fname)
423  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
424  'standard_name', 'projection_x_coordinate'), &
425  this%nc_fname)
426  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
427  'long_name', 'Easting'), this%nc_fname)
428  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, 'bounds', &
429  'mesh_face_xbnds'), this%nc_fname)
430  if (this%wkt /= '') then
431  ! associate with projection
432  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_x, &
433  'grid_mapping', this%gridmap_name), &
434  this%nc_fname)
435  end if
436 
437  ! create mesh x cell bounds variable
438  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_xbnds', nf90_double, &
439  (/this%dim_ids%max_nmesh_face_nodes, &
440  this%dim_ids%nmesh_face/), &
441  this%var_ids%mesh_face_xbnds), &
442  this%nc_fname)
443 
444  ! create mesh y face (cell vertex) variable
445  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_y', nf90_double, &
446  (/this%dim_ids%nmesh_face/), &
447  this%var_ids%mesh_face_y), this%nc_fname)
448 
449  ! assign mesh y face variable attributes
450  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
451  'units', 'm'), this%nc_fname)
452  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
453  'standard_name', 'projection_y_coordinate'), &
454  this%nc_fname)
455  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
456  'long_name', 'Northing'), this%nc_fname)
457  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, 'bounds', &
458  'mesh_face_ybnds'), this%nc_fname)
459 
460  if (this%wkt /= '') then
461  ! associate with projection
462  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_y, &
463  'grid_mapping', this%gridmap_name), &
464  this%nc_fname)
465  end if
466 
467  ! create mesh y cell bounds variable
468  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_ybnds', nf90_double, &
469  (/this%dim_ids%max_nmesh_face_nodes, &
470  this%dim_ids%nmesh_face/), &
471  this%var_ids%mesh_face_ybnds), &
472  this%nc_fname)
473 
474  ! create mesh face nodes variable
475  call nf_verify(nf90_def_var(this%ncid, 'mesh_face_nodes', nf90_int, &
476  (/this%dim_ids%max_nmesh_face_nodes, &
477  this%dim_ids%nmesh_face/), &
478  this%var_ids%mesh_face_nodes), &
479  this%nc_fname)
480 
481  ! assign variable attributes
482  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
483  'cf_role', 'face_node_connectivity'), &
484  this%nc_fname)
485  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
486  'long_name', &
487  'Vertices bounding cell (counterclockwise)'), &
488  this%nc_fname)
489  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
490  '_FillValue', (/nf90_fill_int/)), &
491  this%nc_fname)
492  call nf_verify(nf90_put_att(this%ncid, this%var_ids%mesh_face_nodes, &
493  'start_index', 1), this%nc_fname)
494  end subroutine create_mesh
495 
496  !> @brief define variable chunking
497  !<
498  subroutine ncvar_chunk(ncid, varid, chunk_face, nc_fname)
499  integer(I4B), intent(in) :: ncid
500  integer(I4B), intent(in) :: varid
501  integer(I4B), intent(in) :: chunk_face
502  character(len=*), intent(in) :: nc_fname
503  if (chunk_face > 0) then
504  call nf_verify(nf90_def_var_chunking(ncid, varid, nf90_chunked, &
505  (/chunk_face/)), nc_fname)
506  end if
507  end subroutine ncvar_chunk
508 
509  !> @brief define variable compression
510  !<
511  subroutine ncvar_deflate(ncid, varid, deflate, shuffle, nc_fname)
512  integer(I4B), intent(in) :: ncid
513  integer(I4B), intent(in) :: varid
514  integer(I4B), intent(in) :: deflate
515  integer(I4B), intent(in) :: shuffle
516  character(len=*), intent(in) :: nc_fname
517  if (deflate >= 0) then
518  call nf_verify(nf90_def_var_deflate(ncid, varid, shuffle=shuffle, &
519  deflate=1, deflate_level=deflate), &
520  nc_fname)
521  end if
522  end subroutine ncvar_deflate
523 
524  !> @brief put variable gridmap attributes
525  !<
526  subroutine ncvar_gridmap(ncid, varid, gridmap_name, nc_fname)
527  integer(I4B), intent(in) :: ncid
528  integer(I4B), intent(in) :: varid
529  character(len=*), intent(in) :: gridmap_name
530  character(len=*), intent(in) :: nc_fname
531  if (gridmap_name /= '') then
532  call nf_verify(nf90_put_att(ncid, varid, 'coordinates', &
533  'mesh_face_x mesh_face_y'), nc_fname)
534  call nf_verify(nf90_put_att(ncid, varid, 'grid_mapping', &
535  gridmap_name), nc_fname)
536  end if
537  end subroutine ncvar_gridmap
538 
539  !> @brief put variable internal attributes
540  !<
541  subroutine ncvar_mf6attr(ncid, varid, layer, iper, iaux, nc_tag, nc_fname)
542  integer(I4B), intent(in) :: ncid
543  integer(I4B), intent(in) :: varid
544  integer(I4B), intent(in) :: layer
545  integer(I4B), intent(in) :: iper
546  integer(I4B), intent(in) :: iaux
547  character(len=*), intent(in) :: nc_tag
548  character(len=*), intent(in) :: nc_fname
549  if (nc_tag /= '') then
550  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_input', &
551  nc_tag), nc_fname)
552  if (layer > 0) then
553  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_layer', &
554  layer), nc_fname)
555  end if
556  if (iper > 0) then
557  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_iper', &
558  iper), nc_fname)
559  end if
560  if (iaux > 0) then
561  call nf_verify(nf90_put_att(ncid, varid, 'modflow6_iaux', &
562  iaux), nc_fname)
563  end if
564  end if
565  end subroutine ncvar_mf6attr
566 
567  !> @brief build netcdf variable name
568  !<
569  function export_varname(varname, layer, iper, iaux) result(vname)
570  use inputoutputmodule, only: lowcase
571  character(len=*), intent(in) :: varname
572  integer(I4B), optional, intent(in) :: layer
573  integer(I4B), optional, intent(in) :: iper
574  integer(I4B), optional, intent(in) :: iaux
575  character(len=LINELENGTH) :: vname
576  vname = ''
577  if (varname /= '') then
578  vname = varname
579  call lowcase(vname)
580  if (present(layer)) then
581  if (layer > 0) then
582  write (vname, '(a,i0)') trim(vname)//'_l', layer
583  end if
584  end if
585  if (present(iper)) then
586  if (iper > 0) then
587  write (vname, '(a,i0)') trim(vname)//'_p', iper
588  end if
589  end if
590  if (present(iaux)) then
591  if (iaux > 0) then
592  write (vname, '(a,i0)') trim(vname)//'a', iaux
593  end if
594  end if
595  end if
596  end function export_varname
597 
598 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_warning(msg, substring)
Store warning message.
Definition: Sim.f90:236
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
character(len=maxcharlen) warnmsg
warning message string
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