MODFLOW 6  version 6.6.0.dev0
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, iout)
 initialize More...
 
subroutine mesh_destroy (this)
 initialize 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, iper, iaux, nc_tag, nc_fname)
 put variable internal attributes More...
 
character(len=linelength) function, public export_varname (varname, layer, iper, iaux)
 build netcdf variable name 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 146 of file MeshNCModel.f90.

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)
Here is the call graph for this function:

◆ add_pkg_data()

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

Definition at line 198 of file MeshNCModel.f90.

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)
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 341 of file MeshNCModel.f90.

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)
Here is the call graph for this function:

◆ define_dependent()

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

Definition at line 256 of file MeshNCModel.f90.

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
Here is the call graph for this function:

◆ define_gridmap()

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

Definition at line 317 of file MeshNCModel.f90.

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
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 171 of file MeshNCModel.f90.

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
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:

◆ export_varname()

character(len=linelength) function, public meshmodelmodule::export_varname ( character(len=*), intent(in)  varname,
integer(i4b), intent(in), optional  layer,
integer(i4b), intent(in), optional  iper,
integer(i4b), intent(in), optional  iaux 
)

Definition at line 563 of file MeshNCModel.f90.

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
subroutine, public lowcase(word)
Convert to lower case.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mesh_destroy()

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

Definition at line 136 of file MeshNCModel.f90.

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)
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)  iout 
)
private

Definition at line 98 of file MeshNCModel.f90.

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)
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 492 of file MeshNCModel.f90.

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
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 505 of file MeshNCModel.f90.

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
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 520 of file MeshNCModel.f90.

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
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)  iper,
integer(i4b), intent(in)  iaux,
character(len=*), intent(in)  nc_tag,
character(len=*), intent(in)  nc_fname 
)

Definition at line 535 of file MeshNCModel.f90.

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
Here is the call graph for this function:
Here is the caller graph for this function: