ROMS
Loading...
Searching...
No Matches
get_wetdry_mod Module Reference

Functions/Subroutines

subroutine, public get_wetdry (ng, tile, model, inirec)
 
subroutine, private get_wetdry_nf90 (ng, tile, model, inirec, lbi, ubi, lbj, ubj)
 
subroutine, private get_wetdry_pio (ng, tile, model, inirec, lbi, ubi, lbj, ubj)
 

Function/Subroutine Documentation

◆ get_wetdry()

subroutine, public get_wetdry_mod::get_wetdry ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) inirec )

Definition at line 44 of file get_wetdry.F.

45!***********************************************************************
46!
47! Imported variable declarations.
48!
49 integer, intent(in) :: ng, tile, model, IniRec
50!
51! Local variable declarations.
52!
53 integer :: LBi, UBi, LBj, UBj
54!
55 character (len=*), parameter :: MyFile = &
56 & __FILE__
57!
58!-----------------------------------------------------------------------
59! Read in wetting and drying mask according to IO type.
60!-----------------------------------------------------------------------
61!
62 lbi=bounds(ng)%LBi(tile)
63 ubi=bounds(ng)%UBi(tile)
64 lbj=bounds(ng)%LBj(tile)
65 ubj=bounds(ng)%UBj(tile)
66!
67 SELECT CASE (ini(ng)%IOtype)
68 CASE (io_nf90)
69 CALL get_wetdry_nf90 (ng, tile, model, inirec, &
70 & lbi, ubi, lbj, ubj)
71
72# if defined PIO_LIB && defined DISTRIBUTE
73 CASE (io_pio)
74 CALL get_wetdry_pio (ng, tile, model, inirec, &
75 & lbi, ubi, lbj, ubj)
76# endif
77 CASE DEFAULT
78 IF (master) WRITE (stdout,10) ini(ng)%IOtype
79 exit_flag=2
80 END SELECT
81 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
82!
83 10 FORMAT (' GET_WETDRY - Illegal input file type, io_type = ',i0, &
84 & /,14x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
85!
86 RETURN

References mod_param::bounds, mod_scalars::exit_flag, strings_mod::founderror(), get_wetdry_nf90(), get_wetdry_pio(), mod_iounits::ini, mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_parallel::master, mod_scalars::noerror, and mod_iounits::stdout.

Referenced by ad_initial(), initial(), rp_initial(), and tl_initial().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_wetdry_nf90()

subroutine, private get_wetdry_mod::get_wetdry_nf90 ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) inirec,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 90 of file get_wetdry.F.

92!***********************************************************************
93!
94 USE mod_netcdf
95!
96! Imported variable declarations.
97!
98 integer, intent(in) :: ng, tile, model, IniRec
99 integer, intent(in) :: LBi, UBi, LBj, UBj
100!
101! Local variable declarations.
102!
103 integer :: gtype, i, status, vindex
104 integer :: Vsize(4)
105!
106 real(dp), parameter :: Fscl = 1.0_r8
107
108 real(r8) :: Fmax, Fmin
109!
110 character (len=256) :: ncname
111
112 character (len=*), parameter :: MyFile = &
113 & __FILE__//", get_wetdry_nf90"
114!
115 sourcefile=myfile
116!
117!-----------------------------------------------------------------------
118! Inquire about the contents of restart NetCDF file: Inquire about
119! the dimensions and variables. Check for consistency.
120!-----------------------------------------------------------------------
121!
122 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
123 ncname=ini(ng)%name
124!
125! Open restart (initial) NetCDF for read/write.
126!
127 IF (ini(ng)%ncid.eq.-1) THEN
128 CALL netcdf_open (ng, model, ncname, 1, ini(ng)%ncid)
129 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
130 WRITE (stdout,10) trim(ncname)
131 RETURN
132 END IF
133 END IF
134!
135! Check grid file dimensions for consitency
136!
137 CALL netcdf_check_dim (ng, model, ncname, &
138 & ncid = ini(ng)%ncid)
139 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
140!
141! Inquire about the variables.
142!
143 CALL netcdf_inq_var (ng, model, ncname, &
144 & ncid = ini(ng)%ncid)
145 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
146!
147!-----------------------------------------------------------------------
148! Check if required variables are available.
149!-----------------------------------------------------------------------
150!
151 IF (.not.find_string(var_name,n_var,trim(vname(1,idpwet)), &
152 & vindex)) THEN
153 IF (master) WRITE (stdout,20) trim(vname(1,idpwet)), &
154 & trim(ncname)
155 exit_flag=2
156 RETURN
157 END IF
158 ini(ng)%Vid(idpwet)=vindex
159!
160 IF (.not.find_string(var_name,n_var,trim(vname(1,idrwet)), &
161 & vindex)) THEN
162 IF (master) WRITE (stdout,20) trim(vname(1,idrwet)), &
163 & trim(ncname)
164 exit_flag=2
165 RETURN
166 END IF
167 ini(ng)%Vid(idrwet)=vindex
168!
169 IF (.not.find_string(var_name,n_var,trim(vname(1,iduwet)), &
170 & vindex)) THEN
171 IF (master) WRITE (stdout,20) trim(vname(1,iduwet)), &
172 & trim(ncname)
173 exit_flag=2
174 RETURN
175 END IF
176 ini(ng)%Vid(iduwet)=vindex
177!
178 IF (.not.find_string(var_name,n_var,trim(vname(1,idvwet)), &
179 & vindex)) THEN
180 IF (master) WRITE (stdout,20) trim(vname(1,idvwet)), &
181 & trim(ncname)
182 exit_flag=2
183 RETURN
184 END IF
185 ini(ng)%Vid(idvwet)=vindex
186!
187!-----------------------------------------------------------------------
188! Read in wet/dry mask from rstart file.
189!-----------------------------------------------------------------------
190!
191! Set Vsize to zero to deativate interpolation of input data to model
192! grid in "nf_fread2d".
193!
194 DO i=1,4
195 vsize(i)=0
196 END DO
197!
198! Read in wet/dry mask at PSI-points.
199!
200 gtype=p2dvar
201 status=nf_fread2d(ng, model, ncname, ini(ng)%ncid, &
202 & vname(1,idpwet), ini(ng)%Vid(idpwet), &
203 & inirec, gtype, vsize, &
204 & lbi, ubi, lbj, ubj, &
205 & fscl, fmin, fmax, &
206 & grid(ng) % pmask, &
207 & grid(ng) % pmask_wet)
208 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
209 IF (master) WRITE (stdout,30) trim(vname(1,idpwet)), &
210 & trim(ncname)
211 exit_flag=2
212 ioerror=status
213 RETURN
214 ELSE
215 IF (master) THEN
216 WRITE (stdout,40) trim(vname(2,idpwet)), ng, trim(ncname), &
217 & fmin, fmax
218 END IF
219 END IF
220!
221! Read in wet/dry mask at RHO-points.
222!
223 gtype=r2dvar
224 status=nf_fread2d(ng, model, ncname, ini(ng)%ncid, &
225 & vname(1,idrwet), ini(ng)%Vid(idrwet), &
226 & inirec, gtype, vsize, &
227 & lbi, ubi, lbj, ubj, &
228 & fscl, fmin, fmax, &
229 & grid(ng) % rmask, &
230 & grid(ng) % rmask_wet)
231 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
232 IF (master) WRITE (stdout,30) trim(vname(1,idrwet)), &
233 & trim(ncname)
234 exit_flag=2
235 ioerror=status
236 RETURN
237 ELSE
238 IF (master) THEN
239 WRITE (stdout,40) trim(vname(2,idrwet)), ng, trim(ncname), &
240 & fmin, fmax
241 END IF
242 END IF
243!
244! Read in wet/dry mask at U-points.
245!
246 gtype=u2dvar
247 status=nf_fread2d(ng, model, ncname, ini(ng)%ncid, &
248 & vname(1,iduwet), ini(ng)%Vid(iduwet), &
249 & inirec, gtype, vsize, &
250 & lbi, ubi, lbj, ubj, &
251 & fscl, fmin, fmax, &
252 & grid(ng) % umask, &
253 & grid(ng) % umask_wet)
254 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
255 IF (master) WRITE (stdout,30) trim(vname(1,iduwet)), &
256 & trim(ncname)
257 exit_flag=2
258 ioerror=status
259 RETURN
260 ELSE
261 IF (master) THEN
262 WRITE (stdout,40) trim(vname(2,iduwet)), ng, trim(ncname), &
263 & fmin, fmax
264 END IF
265 END IF
266!
267! Read in wet/dry mask at V-points.
268!
269 gtype=v2dvar
270 status=nf_fread2d(ng, model, ncname, ini(ng)%ncid, &
271 & vname(1,idvwet), ini(ng)%Vid(idvwet), &
272 & inirec, gtype, vsize, &
273 & lbi, ubi, lbj, ubj, &
274 & fscl, fmin, fmax, &
275 & grid(ng) % vmask, &
276 & grid(ng) % vmask_wet)
277 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
278 IF (master) WRITE (stdout,30) trim(vname(1,idvwet)), &
279 & trim(ncname)
280 exit_flag=2
281 ioerror=status
282 RETURN
283 ELSE
284 IF (master) THEN
285 WRITE (stdout,40) trim(vname(2,idvwet)), ng, trim(ncname), &
286 & fmin, fmax
287 END IF
288 END IF
289!
290! Exchange boundary data.
291!
292 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
293 CALL exchange_p2d_tile (ng, tile, &
294 & lbi, ubi, lbj, ubj, &
295 & grid(ng) % pmask_wet)
296 CALL exchange_r2d_tile (ng, tile, &
297 & lbi, ubi, lbj, ubj, &
298 & grid(ng) % rmask_wet)
299 CALL exchange_u2d_tile (ng, tile, &
300 & lbi, ubi, lbj, ubj, &
301 & grid(ng) % umask_wet)
302 CALL exchange_v2d_tile (ng, tile, &
303 & lbi, ubi, lbj, ubj, &
304 & grid(ng) % vmask_wet)
305 END IF
306
307# ifdef DISTRIBUTE
308!
309 CALL mp_exchange2d (ng, tile, model, 4, &
310 & lbi, ubi, lbj, ubj, &
311 & nghostpoints, &
312 & ewperiodic(ng), nsperiodic(ng), &
313 & grid(ng) % pmask_wet, &
314 & grid(ng) % rmask_wet, &
315 & grid(ng) % umask_wet, &
316 & grid(ng) % vmask_wet)
317# endif
318!
319 10 FORMAT (/,' GET_WETDRY_NF90 - unable to open grid NetCDF', &
320 & ' file: ',a)
321 20 FORMAT (/,' GET_WETDRY_NF90 - unable to find grid variable: ',a, &
322 & /,19x,'in NetCDF file: ',a)
323 30 FORMAT (/,' GET_WETDRY_NF90 - error while reading variable: ',a, &
324 & /,19x,'in NetCDF file: ',a)
325 40 FORMAT (2x,'GET_WETDRY_NF90 - ',a,/,23x, &
326 & '(Grid = ',i2.2,', File: ',a,')',/,23x, &
327 & '(Min = ', 1p,e15.8,0p,' Max = ',1p,e15.8,0p,')')
328!
329 RETURN
subroutine, public netcdf_check_dim(ng, model, ncname, ncid)
Definition mod_netcdf.F:538
subroutine, public netcdf_open(ng, model, ncname, omode, ncid)
character(len=100), dimension(mvars) var_name
Definition mod_netcdf.F:169
integer n_var
Definition mod_netcdf.F:152
subroutine, public netcdf_inq_var(ng, model, ncname, ncid, myvarname, searchvar, varid, nvardim, nvaratt)

References mod_scalars::ewperiodic, exchange_2d_mod::exchange_p2d_tile(), exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::exit_flag, strings_mod::find_string(), strings_mod::founderror(), mod_grid::grid, mod_ncparam::idpwet, mod_ncparam::idrwet, mod_ncparam::iduwet, mod_ncparam::idvwet, mod_iounits::ini, mod_iounits::ioerror, mod_parallel::master, mp_exchange_mod::mp_exchange2d(), mod_netcdf::n_var, mod_netcdf::netcdf_check_dim(), mod_netcdf::netcdf_inq_var(), mod_netcdf::netcdf_open(), mod_param::nghostpoints, mod_scalars::noerror, mod_scalars::nsperiodic, mod_param::p2dvar, mod_param::r2dvar, mod_iounits::sourcefile, mod_iounits::stdout, mod_param::u2dvar, mod_param::v2dvar, mod_netcdf::var_name, and mod_ncparam::vname.

Referenced by get_wetdry().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_wetdry_pio()

subroutine, private get_wetdry_mod::get_wetdry_pio ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) inirec,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 335 of file get_wetdry.F.

337!***********************************************************************
338!
340!
341! Imported variable declarations.
342!
343 integer, intent(in) :: ng, tile, model, IniRec
344 integer, intent(in) :: LBi, UBi, LBj, UBj
345!
346! Local variable declarations.
347!
348 integer :: i, status
349 integer :: index_pmask, index_rmask, index_umask, index_vmask
350 integer :: Vsize(4)
351!
352 real(dp), parameter :: Fscl = 1.0_r8
353
354 real(r8) :: Fmax, Fmin
355!
356 character (len=256) :: ncname
357
358 character (len=*), parameter :: MyFile = &
359 & __FILE__//", get_wetdry_pio"
360!
361 TYPE (IO_Desc_t), pointer :: ioDesc
362 TYPE (My_VarDesc), pointer :: pioVar
363!
364 sourcefile=myfile
365!
366!-----------------------------------------------------------------------
367! Inquire about the contents of restart NetCDF file: Inquire about
368! the dimensions and variables. Check for consistency.
369!-----------------------------------------------------------------------
370!
371 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
372 ncname=ini(ng)%name
373!
374! Open restart (initial) NetCDF for read/write.
375!
376 IF (ini(ng)%pioFile%fh.eq.-1) THEN
377 CALL pio_netcdf_open (ng, model, ncname, 1, ini(ng)%pioFile)
378 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
379 WRITE (stdout,10) trim(ncname)
380 RETURN
381 END IF
382 END IF
383!
384! Check grid file dimensions for consitency
385!
386 CALL pio_netcdf_check_dim (ng, model, ncname, &
387 & piofile = ini(ng)%pioFile)
388 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
389!
390! Inquire about the variables.
391!
392 CALL netcdf_inq_var (ng, model, ncname, &
393 & piofile = ini(ng)%pioFile)
394 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
395!
396!-----------------------------------------------------------------------
397! Check if required variables are available.
398!-----------------------------------------------------------------------
399!
400 IF (.not.find_string(var_name,n_var,trim(vname(1,idpwet)), &
401 & index_pmask)) THEN
402 IF (master) WRITE (stdout,20) trim(vname(1,idpwet)), &
403 & trim(ncname)
404 exit_flag=2
405 RETURN
406 END IF
407!
408 IF (.not.find_string(var_name,n_var,trim(vname(1,idrwet)), &
409 & index_rmask)) THEN
410 IF (master) WRITE (stdout,20) trim(vname(1,idrwet)), &
411 & trim(ncname)
412 exit_flag=2
413 RETURN
414 END IF
415!
416 IF (.not.find_string(var_name,n_var,trim(vname(1,iduwet)), &
417 & index_umask)) THEN
418 IF (master) WRITE (stdout,20) trim(vname(1,iduwet)), &
419 & trim(ncname)
420 exit_flag=2
421 RETURN
422 END IF
423!
424 IF (.not.find_string(var_name,n_var,trim(vname(1,idvwet)), &
425 & index_vmask)) THEN
426 IF (master) WRITE (stdout,20) trim(vname(1,idvwet)), &
427 & trim(ncname)
428 exit_flag=2
429 RETURN
430 END IF
431!
432!-----------------------------------------------------------------------
433! Read in wet/dry mask from rstart file.
434!-----------------------------------------------------------------------
435!
436! Set Vsize to zero to deativate interpolation of input data to model
437! grid in "nf_fread2d".
438!
439 DO i=1,4
440 vsize(i)=0
441 END DO
442!
443! Read in wet/dry mask at PSI-points.
444!
445 piovar%vd => var_desc(index_pmask)
446 IF (kind(grid(ng)%pmask_wet).eq.8) THEN
447 piovar%dkind=pio_double
448 iodesc => iodesc_dp_p2dvar(ng)
449 ELSE
450 piovar%dkind=pio_real
451 iodesc => iodesc_sp_p2dvar(ng)
452 END IF
453 piovar%gtype=p2dvar
454!
455 status=nf_fread2d(ng, model, ncname, ini(ng)%pioFile, &
456 & vname(1,idpwet), piovar, &
457 & inirec, iodesc, vsize, &
458 & lbi, ubi, lbj, ubj, &
459 & fscl, fmin, fmax, &
460 & grid(ng) % pmask, &
461 & grid(ng) % pmask_wet)
462 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
463 IF (master) WRITE (stdout,30) trim(vname(1,idpwet)), &
464 & trim(ncname)
465 exit_flag=2
466 ioerror=status
467 RETURN
468 ELSE
469 IF (master) THEN
470 WRITE (stdout,40) trim(vname(2,idpwet)), ng, trim(ncname), &
471 & fmin, fmax
472 END IF
473 END IF
474!
475! Read in wet/dry mask at RHO-points.
476!
477 piovar%vd => var_desc(index_rmask)
478 IF (kind(grid(ng)%rmask_wet).eq.8) THEN
479 piovar%dkind=pio_double
480 iodesc => iodesc_dp_r2dvar(ng)
481 ELSE
482 piovar%dkind=pio_real
483 iodesc => iodesc_sp_r2dvar(ng)
484 END IF
485 piovar%gtype=r2dvar
486!
487 status=nf_fread2d(ng, model, ncname, ini(ng)%pioFile, &
488 & vname(1,idrwet), piovar, &
489 & inirec, iodesc, vsize, &
490 & lbi, ubi, lbj, ubj, &
491 & fscl, fmin, fmax, &
492 & grid(ng) % rmask, &
493 & grid(ng) % rmask_wet)
494 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
495 IF (master) WRITE (stdout,30) trim(vname(1,idrwet)), &
496 & trim(ncname)
497 exit_flag=2
498 ioerror=status
499 RETURN
500 ELSE
501 IF (master) THEN
502 WRITE (stdout,40) trim(vname(2,idrwet)), ng, trim(ncname), &
503 & fmin, fmax
504 END IF
505 END IF
506!
507! Read in wet/dry mask at U-points.
508!
509 piovar%vd => var_desc(index_umask)
510 IF (kind(grid(ng)%umask_wet).eq.8) THEN
511 piovar%dkind=pio_double
512 iodesc => iodesc_dp_u2dvar(ng)
513 ELSE
514 piovar%dkind=pio_real
515 iodesc => iodesc_sp_u2dvar(ng)
516 END IF
517 piovar%gtype=u2dvar
518!
519 status=nf_fread2d(ng, model, ncname, ini(ng)%pioFile, &
520 & vname(1,iduwet), piovar, &
521 & inirec, iodesc, vsize, &
522 & lbi, ubi, lbj, ubj, &
523 & fscl, fmin, fmax, &
524 & grid(ng) % umask, &
525 & grid(ng) % umask_wet)
526 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
527 IF (master) WRITE (stdout,30) trim(vname(1,iduwet)), &
528 & trim(ncname)
529 exit_flag=2
530 ioerror=status
531 RETURN
532 ELSE
533 IF (master) THEN
534 WRITE (stdout,40) trim(vname(2,iduwet)), ng, trim(ncname), &
535 & fmin, fmax
536 END IF
537 END IF
538!
539! Read in wet/dry mask at V-points.
540!
541 piovar%vd => var_desc(index_vmask)
542 IF (kind(grid(ng)%vmask_wet).eq.8) THEN
543 piovar%dkind=pio_double
544 iodesc => iodesc_dp_v2dvar(ng)
545 ELSE
546 piovar%dkind=pio_real
547 iodesc => iodesc_sp_v2dvar(ng)
548 END IF
549 piovar%gtype=v2dvar
550!
551 status=nf_fread2d(ng, model, ncname, ini(ng)%pioFile, &
552 & vname(1,idvwet), piovar, &
553 & inirec, iodesc, vsize, &
554 & lbi, ubi, lbj, ubj, &
555 & fscl, fmin, fmax, &
556 & grid(ng) % vmask, &
557 & grid(ng) % vmask_wet)
558 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
559 IF (master) WRITE (stdout,30) trim(vname(1,idvwet)), &
560 & trim(ncname)
561 exit_flag=2
562 ioerror=status
563 RETURN
564 ELSE
565 IF (master) THEN
566 WRITE (stdout,40) trim(vname(2,idvwet)), ng, trim(ncname), &
567 & fmin, fmax
568 END IF
569 END IF
570!
571! Exchange boundary data.
572!
573 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
574 CALL exchange_p2d_tile (ng, tile, &
575 & lbi, ubi, lbj, ubj, &
576 & grid(ng) % pmask_wet)
577 CALL exchange_r2d_tile (ng, tile, &
578 & lbi, ubi, lbj, ubj, &
579 & grid(ng) % rmask_wet)
580 CALL exchange_u2d_tile (ng, tile, &
581 & lbi, ubi, lbj, ubj, &
582 & grid(ng) % umask_wet)
583 CALL exchange_v2d_tile (ng, tile, &
584 & lbi, ubi, lbj, ubj, &
585 & grid(ng) % vmask_wet)
586 END IF
587!
588 CALL mp_exchange2d (ng, tile, model, 4, &
589 & lbi, ubi, lbj, ubj, &
590 & nghostpoints, &
591 & ewperiodic(ng), nsperiodic(ng), &
592 & grid(ng) % pmask_wet, &
593 & grid(ng) % rmask_wet, &
594 & grid(ng) % umask_wet, &
595 & grid(ng) % vmask_wet)
596!
597 10 FORMAT (/,' GET_WETDRY_PIO - unable to open grid NetCDF', &
598 & ' file: ',a)
599 20 FORMAT (/,' GET_WETDRY_PIO - unable to find grid variable: ',a, &
600 & /,18x,'in NetCDF file: ',a)
601 30 FORMAT (/,' GET_WETDRY_PIO - error while reading variable: ',a, &
602 & /,18x,'in NetCDF file: ',a)
603 40 FORMAT (2x,'GET_WETDRY_PIO - ',a,/,23x, &
604 & '(Grid = ',i2.2,', File: ',a,')',/,23x, &
605 & '(Min = ', 1p,e15.8,0p,' Max = ',1p,e15.8,0p,')')
606!
607 RETURN
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dvar
type(var_desc_t), dimension(:), pointer var_desc
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dvar
subroutine, public pio_netcdf_open(ng, model, ncname, omode, piofile)
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_p2dvar
subroutine, public pio_netcdf_check_dim(ng, model, ncname, piofile)
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_p2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dvar

References mod_scalars::ewperiodic, exchange_2d_mod::exchange_p2d_tile(), exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::exit_flag, strings_mod::find_string(), strings_mod::founderror(), mod_grid::grid, mod_ncparam::idpwet, mod_ncparam::idrwet, mod_ncparam::iduwet, mod_ncparam::idvwet, mod_iounits::ini, mod_pio_netcdf::iodesc_dp_p2dvar, mod_pio_netcdf::iodesc_dp_r2dvar, mod_pio_netcdf::iodesc_dp_u2dvar, mod_pio_netcdf::iodesc_dp_v2dvar, mod_pio_netcdf::iodesc_sp_p2dvar, mod_pio_netcdf::iodesc_sp_r2dvar, mod_pio_netcdf::iodesc_sp_u2dvar, mod_pio_netcdf::iodesc_sp_v2dvar, mod_iounits::ioerror, mod_parallel::master, mp_exchange_mod::mp_exchange2d(), mod_param::nghostpoints, mod_scalars::noerror, mod_scalars::nsperiodic, mod_param::p2dvar, mod_pio_netcdf::pio_netcdf_check_dim(), mod_pio_netcdf::pio_netcdf_open(), mod_param::r2dvar, mod_iounits::sourcefile, mod_iounits::stdout, mod_param::u2dvar, mod_param::v2dvar, mod_pio_netcdf::var_desc, and mod_ncparam::vname.

Referenced by get_wetdry().

Here is the call graph for this function:
Here is the caller graph for this function: