MODFLOW 6  version 6.7.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 152 of file MeshNCModel.f90.

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

◆ add_pkg_data()

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

Definition at line 204 of file MeshNCModel.f90.

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

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

◆ define_dependent()

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

Definition at line 262 of file MeshNCModel.f90.

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

◆ define_gridmap()

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

Definition at line 323 of file MeshNCModel.f90.

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

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

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

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

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

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

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

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