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

Functions/Subroutines

subroutine, public state_read (ng, tile, model, iotype, lbi, ubi, lbj, ubj, lbij, ubij, lout, rec, nopen, ncid, piofile, ncname, rmask, umask, vmask, s_t_obc, s_u_obc, s_v_obc, s_ubar_obc, s_vbar_obc, s_zeta_obc, s_ustr, s_vstr, s_tflux, s_t, s_u, s_v, s_zeta)
 
subroutine, private state_read_nf90 (ng, tile, model, lreport, lbi, ubi, lbj, ubj, lbij, ubij, lout, rec, nopen, ncid, ncname, rmask, umask, vmask, s_t_obc, s_u_obc, s_v_obc, s_ubar_obc, s_vbar_obc, s_zeta_obc, s_ustr, s_vstr, s_tflux, s_t, s_u, s_v, s_zeta)
 
subroutine, private state_read_pio (ng, tile, model, lreport, lbi, ubi, lbj, ubj, lbij, ubij, lout, rec, nopen, piofile, ncname, ifdef masking
 

Function/Subroutine Documentation

◆ state_read()

subroutine, public state_read_mod::state_read ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) iotype,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lout,
integer, intent(in) rec,
integer, intent(in) nopen,
integer, intent(inout) ncid,
type (file_desc_t), intent(inout) piofile,
character (len=*), intent(in) ncname,
real(r8), dimension(lbi:,lbj:), intent(in) rmask,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbij:,:,:,:,:,:), intent(inout) s_t_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s_v_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_zeta_obc,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_ustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_vstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s_t,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_v,
real(r8), dimension(lbi:,lbj:,:), intent(inout) s_zeta )

Definition at line 107 of file state_read.F.

137!***********************************************************************
138!
139! Imported variable declarations.
140!
141 integer, intent(in) :: ng, tile, model, IOtype
142 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
143 integer, intent(in) :: Lout, rec, nopen
144
145 integer, intent(inout) :: ncid
146
147#if defined PIO_LIB && defined DISTRIBUTE
148!
149 TYPE (File_desc_t), intent(inout) :: pioFile
150#endif
151!
152 character (len=*), intent(in) :: ncname
153!
154#ifdef ASSUMED_SHAPE
155# ifdef MASKING
156 real(r8), intent(in) :: rmask(LBi:,LBj:)
157 real(r8), intent(in) :: umask(LBi:,LBj:)
158 real(r8), intent(in) :: vmask(LBi:,LBj:)
159# endif
160# ifdef ADJUST_BOUNDARY
161# ifdef SOLVE3D
162 real(r8), intent(inout) :: s_t_obc(LBij:,:,:,:,:,:)
163 real(r8), intent(inout) :: s_u_obc(LBij:,:,:,:,:)
164 real(r8), intent(inout) :: s_v_obc(LBij:,:,:,:,:)
165# endif
166 real(r8), intent(inout) :: s_ubar_obc(LBij:,:,:,:)
167 real(r8), intent(inout) :: s_vbar_obc(LBij:,:,:,:)
168 real(r8), intent(inout) :: s_zeta_obc(LBij:,:,:,:)
169# endif
170# ifdef ADJUST_WSTRESS
171 real(r8), intent(inout) :: s_ustr(LBi:,LBj:,:,:)
172 real(r8), intent(inout) :: s_vstr(LBi:,LBj:,:,:)
173# endif
174# ifdef SOLVE3D
175# ifdef ADJUST_STFLUX
176 real(r8), intent(inout) :: s_tflux(LBi:,LBj:,:,:,:)
177# endif
178 real(r8), intent(inout) :: s_t(LBi:,LBj:,:,:,:)
179 real(r8), intent(inout) :: s_u(LBi:,LBj:,:,:)
180 real(r8), intent(inout) :: s_v(LBi:,LBj:,:,:)
181# else
182 real(r8), intent(inout) :: s_ubar(LBi:,LBj:,:)
183 real(r8), intent(inout) :: s_vbar(LBi:,LBj:,:)
184# endif
185 real(r8), intent(inout) :: s_zeta(LBi:,LBj:,:)
186
187#else
188
189# ifdef MASKING
190 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
191 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
192 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
193# endif
194# ifdef ADJUST_BOUNDARY
195# ifdef SOLVE3D
196 real(r8), intent(inout) :: s_t_obc(LBij:UBij,N(ng),4, &
197 & Nbrec(ng),2,NT(ng))
198 real(r8), intent(inout) :: s_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
199 real(r8), intent(inout) :: s_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
200# endif
201 real(r8), intent(inout) :: s_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
202 real(r8), intent(inout) :: s_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
203 real(r8), intent(inout) :: s_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
204# endif
205# ifdef ADJUST_WSTRESS
206 real(r8), intent(inout) :: s_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
207 real(r8), intent(inout) :: s_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
208# endif
209# ifdef SOLVE3D
210# ifdef ADJUST_STFLUX
211 real(r8), intent(inout) :: s_tflux(LBi:UBi,LBj:UBj, &
212 & Nfrec(ng),2,NT(ng))
213# endif
214 real(r8), intent(inout) :: s_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
215 real(r8), intent(inout) :: s_u(LBi:UBi,LBj:UBj,N(ng),2)
216 real(r8), intent(inout) :: s_v(LBi:UBi,LBj:UBj,N(ng),2)
217# else
218 real(r8), intent(inout) :: s_ubar(LBi:UBi,LBj:UBj,:)
219 real(r8), intent(inout) :: s_vbar(LBi:UBi,LBj:UBj,:)
220# endif
221 real(r8), intent(inout) :: s_zeta(LBi:UBi,LBj:UBj,:)
222#endif
223!
224! Local variable declarations.
225!
226 logical :: Lreport = .false.
227!
228 character (len=*), parameter :: MyFile = &
229 & __FILE__
230!
231 sourcefile=myfile
232!
233!-----------------------------------------------------------------------
234! Read in requested state vector
235!-----------------------------------------------------------------------
236!
237 SELECT CASE (iotype)
238 CASE (io_nf90)
239 CALL state_read_nf90 (ng, tile, model, lreport, &
240 & lbi, ubi, lbj, ubj, lbij, ubij, &
241 & lout, rec, &
242 & nopen, ncid, ncname, &
243#ifdef MASKING
244 & rmask, umask, vmask, &
245#endif
246#ifdef ADJUST_BOUNDARY
247# ifdef SOLVE3D
248 & s_t_obc, s_u_obc, s_v_obc, &
249# endif
250 & s_ubar_obc, s_vbar_obc, &
251 & s_zeta_obc, &
252#endif
253#ifdef ADJUST_WSTRESS
254 & s_ustr, s_vstr, &
255#endif
256#ifdef SOLVE3D
257# ifdef ADJUST_STFLUX
258 & s_tflux, &
259# endif
260 & s_t, s_u, s_v, &
261#else
262 & s_ubar, s_vbar, &
263#endif
264 & s_zeta)
265
266#if defined PIO_LIB && defined DISTRIBUTE
267 CASE (io_pio)
268 CALL state_read_pio (ng, tile, model, lreport, &
269 & lbi, ubi, lbj, ubj, lbij, ubij, &
270 & lout, rec, &
271 & nopen, piofile, ncname, &
272# ifdef MASKING
273 & rmask, umask, vmask, &
274# endif
275# ifdef ADJUST_BOUNDARY
276# ifdef SOLVE3D
277 & s_t_obc, s_u_obc, s_v_obc, &
278# endif
279 & s_ubar_obc, s_vbar_obc, &
280 & s_zeta_obc, &
281# endif
282# ifdef ADJUST_WSTRESS
283 & s_ustr, s_vstr, &
284# endif
285# ifdef SOLVE3D
286# ifdef ADJUST_STFLUX
287 & s_tflux, &
288# endif
289 & s_t, s_u, s_v, &
290# else
291 & s_ubar, s_vbar, &
292# endif
293 & s_zeta)
294#endif
295 CASE DEFAULT
296 IF (master) WRITE (stdout,10) iotype
297 exit_flag=2
298 END SELECT
299 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
300!
301 10 FORMAT (' STATE_READ - Illegal inpput file type, io_type = ',i0, &
302 & /,14x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
303!
304 RETURN

References mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_parallel::master, mod_scalars::noerror, mod_iounits::sourcefile, state_read_nf90(), state_read_pio(), and mod_iounits::stdout.

Referenced by posterior_mod::analysis_error(), cgradient_mod::hessian(), cgradient_mod::hessian_evecs(), ini_lanczos_mod::ini_lanczos_tile(), cgradient_mod::lanczos(), posterior_mod::lanczos(), cgradient_mod::new_cost(), cgradient_mod::new_gradient(), posterior_mod::posterior_eofs(), posterior_var_mod::posterior_var_tile(), cgradient_mod::precond(), inner2state_mod::tl_inner2state_tile(), and cgradient_mod::tl_new_state().

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

◆ state_read_nf90()

subroutine, private state_read_mod::state_read_nf90 ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
logical, intent(in) lreport,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lout,
integer, intent(in) rec,
integer, intent(in) nopen,
integer, intent(inout) ncid,
character (len=*), intent(in) ncname,
real(r8), dimension(lbi:,lbj:), intent(in) rmask,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbij:,:,:,:,:,:), intent(inout) s_t_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s_u_obc,
real(r8), dimension(lbij:,:,:,:,:), intent(inout) s_v_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_ubar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_vbar_obc,
real(r8), dimension(lbij:,:,:,:), intent(inout) s_zeta_obc,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_ustr,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_vstr,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s_tflux,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) s_t,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) s_v,
real(r8), dimension(lbi:,lbj:,:), intent(inout) s_zeta )
private

Definition at line 308 of file state_read.F.

334!***********************************************************************
335!
336 USE mod_netcdf
337!
338! Imported variable declarations.
339!
340 logical, intent(in) :: Lreport
341!
342 integer, intent(in) :: ng, tile, model
343 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
344 integer, intent(in) :: Lout, rec, nopen
345
346 integer, intent(inout) :: ncid
347!
348 character (len=*), intent(in) :: ncname
349!
350#ifdef ASSUMED_SHAPE
351# ifdef MASKING
352 real(r8), intent(in) :: rmask(LBi:,LBj:)
353 real(r8), intent(in) :: umask(LBi:,LBj:)
354 real(r8), intent(in) :: vmask(LBi:,LBj:)
355# endif
356# ifdef ADJUST_BOUNDARY
357# ifdef SOLVE3D
358 real(r8), intent(inout) :: s_t_obc(LBij:,:,:,:,:,:)
359 real(r8), intent(inout) :: s_u_obc(LBij:,:,:,:,:)
360 real(r8), intent(inout) :: s_v_obc(LBij:,:,:,:,:)
361# endif
362 real(r8), intent(inout) :: s_ubar_obc(LBij:,:,:,:)
363 real(r8), intent(inout) :: s_vbar_obc(LBij:,:,:,:)
364 real(r8), intent(inout) :: s_zeta_obc(LBij:,:,:,:)
365# endif
366# ifdef ADJUST_WSTRESS
367 real(r8), intent(inout) :: s_ustr(LBi:,LBj:,:,:)
368 real(r8), intent(inout) :: s_vstr(LBi:,LBj:,:,:)
369# endif
370# ifdef SOLVE3D
371# ifdef ADJUST_STFLUX
372 real(r8), intent(inout) :: s_tflux(LBi:,LBj:,:,:,:)
373# endif
374 real(r8), intent(inout) :: s_t(LBi:,LBj:,:,:,:)
375 real(r8), intent(inout) :: s_u(LBi:,LBj:,:,:)
376 real(r8), intent(inout) :: s_v(LBi:,LBj:,:,:)
377# else
378 real(r8), intent(inout) :: s_ubar(LBi:,LBj:,:)
379 real(r8), intent(inout) :: s_vbar(LBi:,LBj:,:)
380# endif
381 real(r8), intent(inout) :: s_zeta(LBi:,LBj:,:)
382
383#else
384
385# ifdef MASKING
386 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
387 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
388 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
389# endif
390# ifdef ADJUST_BOUNDARY
391# ifdef SOLVE3D
392 real(r8), intent(inout) :: s_t_obc(LBij:UBij,N(ng),4, &
393 & Nbrec(ng),2,NT(ng))
394 real(r8), intent(inout) :: s_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
395 real(r8), intent(inout) :: s_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
396# endif
397 real(r8), intent(inout) :: s_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
398 real(r8), intent(inout) :: s_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
399 real(r8), intent(inout) :: s_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
400# endif
401# ifdef ADJUST_WSTRESS
402 real(r8), intent(inout) :: s_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
403 real(r8), intent(inout) :: s_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
404# endif
405# ifdef SOLVE3D
406# ifdef ADJUST_STFLUX
407 real(r8), intent(inout) :: s_tflux(LBi:UBi,LBj:UBj, &
408 & Nfrec(ng),2,NT(ng))
409# endif
410 real(r8), intent(inout) :: s_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
411 real(r8), intent(inout) :: s_u(LBi:UBi,LBj:UBj,N(ng),2)
412 real(r8), intent(inout) :: s_v(LBi:UBi,LBj:UBj,N(ng),2)
413# else
414 real(r8), intent(inout) :: s_ubar(LBi:UBi,LBj:UBj,:)
415 real(r8), intent(inout) :: s_vbar(LBi:UBi,LBj:UBj,:)
416# endif
417 real(r8), intent(inout) :: s_zeta(LBi:UBi,LBj:UBj,:)
418#endif
419!
420! Local variable declarations.
421!
422 integer :: Sstr, Send
423 integer :: i, j, k
424 integer :: ifield, itrc
425 integer :: gtype, my_ncid, status, varid
426
427 integer, dimension(4) :: Vsize
428!
429 real(r8) :: Fmin, Fmax
430 real(dp) :: stime, scale
431!
432 character (len=15) :: Tstring
433 character (len=22) :: t_code
434
435 character (len=*), parameter :: MyFile = &
436 & __FILE__//", state_read_nf90"
437
438#include "set_bounds.h"
439!
440 sourcefile=myfile
441!
442!-----------------------------------------------------------------------
443! Read in requested model state record. Load data into state array
444! index Lout.
445!-----------------------------------------------------------------------
446#ifdef PROFILE
447!
448! Turn on time wall clock.
449!
450 CALL wclock_on (ng, model, 80, __line__, myfile)
451#endif
452!
453! Determine file and variables ids.
454!
455 IF ((nopen.gt.0).or.(ncid.eq.-1)) THEN
456 CALL netcdf_open (ng, model, ncname, 0, my_ncid)
457 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
458 ncid=my_ncid
459 ELSE
460 my_ncid=ncid
461 END IF
462!
463 DO i=1,4
464 vsize(i)=0
465 END DO
466
467#ifdef SP4DVAR
468!
469! Synchronize NetCDF data to disk. It allows processes open for reading
470! data to update the number of records written by other processes.
471! The writer flushes buffers to disk, and the reader makes sure that
472! it is reading from disk rather than from previously cached buffers.
473!
474 CALL netcdf_sync (ng, model, ncname, my_ncid)
475 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
476#endif
477!
478! Read in time.
479!
480 CALL netcdf_get_time (ng, model, ncname, vname(1,idtime), &
481 & rclock%DateNumber, stime, &
482 & my_ncid, (/rec/), (/1/))
483 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
484!
485! Report information.
486!
487 IF (master) THEN
488 CALL time_string (stime, t_code)
489 sstr=scan(calledfrom,'/',back=.true.)+1
490 send=len_trim(calledfrom)
491 WRITE (tstring,'(f15.4)') stime*sec2day
492 WRITE (stdout,10) 'Reading state fields,', t_code, &
493 & ng, tstring, trim(ncname), rec, lout, &
494 & calledfrom(sstr:send)
495 END IF
496!
497! Read in free-surface.
498!
499 gtype=r2dvar
500 scale=1.0_dp
501 CALL netcdf_inq_varid (ng, model, ncname, vname(1,idfsur), &
502 & my_ncid, varid)
503 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
504
505 status=nf_fread2d(ng, model, ncname, my_ncid, &
506 & vname(1,idfsur), varid, &
507 & rec, gtype, vsize, &
508 & lbi, ubi, lbj, ubj, &
509 & scale, fmin, fmax, &
510#ifdef MASKING
511 & rmask, &
512#endif
513 & s_zeta(:,:,lout))
514 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
515 IF (master) THEN
516 WRITE (stdout,30) trim(vname(1,idfsur)), rec, trim(ncname)
517 END IF
518 exit_flag=3
519 ioerror=status
520 RETURN
521 ELSE
522 IF (master.and.lreport) THEN
523 WRITE (stdout,40) trim(vname(2,idfsur)), fmin, fmax
524 END IF
525 END IF
526
527#ifdef ADJUST_BOUNDARY
528!
529! Read in free-surface open boundaries.
530!
531 IF (any(lobc(:,isfsur,ng))) THEN
532 ifield=idsbry(isfsur)
533 gtype=r2dvar
534 scale=1.0_dp
535 CALL netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
536 & my_ncid, varid)
537 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
538
539 status=nf_fread2d_bry(ng, model, ncname, my_ncid, &
540 & vname(1,ifield), varid, &
541 & rec, gtype, &
542 & lbij, ubij, nbrec(ng), &
543 & scale, fmin, fmax, &
544 & s_zeta_obc(:,:,:,lout))
545 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
546 IF (master) THEN
547 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
548 END IF
549 exit_flag=3
550 ioerror=status
551 RETURN
552 ELSE
553 IF (master.and.lreport) THEN
554 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
555 END IF
556 END IF
557 END IF
558#endif
559
560#ifndef SOLVE3D
561!
562! Read in 2D U-momentum component.
563!
564 gtype=u2dvar
565 scale=1.0_dp
566 CALL netcdf_inq_varid (ng, model, ncname, vname(1,idubar), &
567 & my_ncid, varid)
568 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
569
570 status=nf_fread2d(ng, model, ncname, my_ncid, &
571 & vname(1,idubar), varid, &
572 & rec, gtype, vsize, &
573 & lbi, ubi, lbj, ubj, &
574 & scale, fmin, fmax, &
575# ifdef MASKING
576 & umask, &
577# endif
578 & s_ubar(:,:,lout))
579 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
580 IF (master) THEN
581 WRITE (stdout,30) trim(vname(1,idubar)), rec, trim(ncname)
582 END IF
583 exit_flag=3
584 ioerror=status
585 RETURN
586 ELSE
587 IF (master.and.lreport) THEN
588 WRITE (stdout,40) trim(vname(2,idubar)), fmin, fmax
589 END IF
590 END IF
591!
592! Read in 2D V-momentum component.
593!
594 gtype=v2dvar
595 scale=1.0_dp
596 CALL netcdf_inq_varid (ng, model, ncname, vname(1,idvbar), &
597 & my_ncid, varid)
598 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
599
600 status=nf_fread2d(ng, model, ncname, my_ncid, &
601 & vname(1,idvbar), varid, &
602 & rec, gtype, vsize, &
603 & lbi, ubi, lbj, ubj, &
604 & scale, fmin, fmax, &
605# ifdef MASKING
606 & vmask, &
607# endif
608 & s_vbar(:,:,lout))
609 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
610 IF (master) THEN
611 WRITE (stdout,30) trim(vname(1,idvbar)), rec, trim(ncname)
612 END IF
613 exit_flag=3
614 ioerror=status
615 RETURN
616 ELSE
617 IF (master.and.lreport) THEN
618 WRITE (stdout,40) trim(vname(2,idvbar)), fmin, fmax
619 END IF
620 END IF
621#endif
622
623#ifdef ADJUST_BOUNDARY
624!
625! Read in 2D U-momentum component open boundaries.
626!
627 IF (any(lobc(:,isubar,ng))) THEN
628 ifield=idsbry(isubar)
629 gtype=u2dvar
630 scale=1.0_dp
631 CALL netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
632 & my_ncid, varid)
633 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
634
635 status=nf_fread2d_bry(ng, model, ncname, my_ncid, &
636 & vname(1,ifield), varid, &
637 & rec, gtype, &
638 & lbij, ubij, nbrec(ng), &
639 & scale, fmin, fmax, &
640 & s_ubar_obc(:,:,:,lout))
641 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
642 IF (master) THEN
643 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
644 END IF
645 exit_flag=3
646 ioerror=status
647 RETURN
648 ELSE
649 IF (master.and.lreport) THEN
650 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
651 END IF
652 END IF
653 END IF
654!
655! Read in 2D V-momentum component open boundaries.
656!
657 IF (any(lobc(:,isvbar,ng))) THEN
658 ifield=idsbry(isvbar)
659 gtype=v2dvar
660 scale=1.0_dp
661 CALL netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
662 & my_ncid, varid)
663 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
664
665 status=nf_fread2d_bry(ng, model, ncname, my_ncid, &
666 & vname(1,ifield), varid, &
667 & rec, gtype, &
668 & lbij, ubij, nbrec(ng), &
669 & scale, fmin, fmax, &
670 & s_vbar_obc(:,:,:,lout))
671 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
672 IF (master) THEN
673 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
674 END IF
675 exit_flag=3
676 ioerror=status
677 RETURN
678 ELSE
679 IF (master.and.lreport) THEN
680 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
681 END IF
682 END IF
683 END IF
684#endif
685
686#ifdef ADJUST_WSTRESS
687!
688! Read surface momentum stress.
689!
690 gtype=u3dvar
691 scale=1.0_dp
692 CALL netcdf_inq_varid (ng, model, ncname, vname(1,idusms), &
693 & my_ncid, varid)
694 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
695
696 status=nf_fread3d(ng, model, ncname, my_ncid, &
697 & vname(1,idusms), varid, &
698 & rec, gtype, vsize, &
699 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
700 & scale, fmin, fmax, &
701# ifdef MASKING
702 & umask, &
703# endif
704 & s_ustr(:,:,:,lout))
705 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
706 IF (master) THEN
707 WRITE (stdout,30) trim(vname(1,idusms)), rec, trim(ncname)
708 END IF
709 exit_flag=3
710 ioerror=status
711 RETURN
712 ELSE
713 IF (master.and.lreport) THEN
714 WRITE (stdout,40) trim(vname(2,idusms)), fmin, fmax
715 END IF
716 END IF
717!
718 gtype=v3dvar
719 scale=1.0_dp
720 CALL netcdf_inq_varid (ng, model, ncname, vname(1,idvsms), &
721 & my_ncid, varid)
722 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
723
724 status=nf_fread3d(ng, model, ncname, my_ncid, &
725 & vname(1,idvsms), varid, &
726 & rec, gtype, vsize, &
727 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
728 & scale, fmin, fmax, &
729# ifdef MASKING
730 & vmask, &
731# endif
732 & s_vstr(:,:,:,lout))
733 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
734 IF (master) THEN
735 WRITE (stdout,30) trim(vname(1,idvsms)), rec, trim(ncname)
736 END IF
737 exit_flag=3
738 ioerror=status
739 RETURN
740 ELSE
741 IF (master.and.lreport) THEN
742 WRITE (stdout,40) trim(vname(2,idvsms)), fmin, fmax
743 END IF
744 END IF
745#endif
746
747#ifdef SOLVE3D
748!
749! Read in 3D U-momentum component.
750!
751 gtype=u3dvar
752 scale=1.0_dp
753 CALL netcdf_inq_varid (ng, model, ncname, vname(1,iduvel), &
754 & my_ncid, varid)
755 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
756
757 status=nf_fread3d(ng, model, ncname, my_ncid, &
758 & vname(1,iduvel), varid, &
759 & rec, gtype, vsize, &
760 & lbi, ubi, lbj, ubj, 1, n(ng), &
761 & scale, fmin, fmax, &
762# ifdef MASKING
763 & umask, &
764# endif
765 & s_u(:,:,:,lout))
766 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
767 IF (master) THEN
768 WRITE (stdout,30) trim(vname(1,iduvel)), rec, trim(ncname)
769 END IF
770 exit_flag=3
771 ioerror=status
772 RETURN
773 ELSE
774 IF (master.and.lreport) THEN
775 WRITE (stdout,40) trim(vname(2,iduvel)), fmin, fmax
776 END IF
777 END IF
778
779# ifdef ADJUST_BOUNDARY
780!
781! Read in 3D U-momentum component open boundaries.
782!
783 IF (any(lobc(:,isuvel,ng))) THEN
784 ifield=idsbry(isuvel)
785 gtype=u3dvar
786 scale=1.0_dp
787 CALL netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
788 & my_ncid, varid)
789 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
790
791 status=nf_fread3d_bry(ng, model, ncname, my_ncid, &
792 & vname(1,ifield), varid, &
793 & rec, gtype, &
794 & lbij, ubij, 1, n(ng), nbrec(ng), &
795 & scale, fmin, fmax, &
796 & s_u_obc(:,:,:,:,lout))
797 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
798 IF (master) THEN
799 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
800 END IF
801 exit_flag=3
802 ioerror=status
803 RETURN
804 ELSE
805 IF (master.and.lreport) THEN
806 WRITE (stdout,40) trim(vname(1,ifield)), fmin, fmax
807 END IF
808 END IF
809 END IF
810# endif
811!
812! Read in 3D V-momentum component.
813!
814 gtype=v3dvar
815 scale=1.0_dp
816 CALL netcdf_inq_varid (ng, model, ncname, vname(1,idvvel), &
817 & my_ncid, varid)
818 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
819
820 status=nf_fread3d(ng, model, ncname, my_ncid, &
821 & vname(1,idvvel), varid, &
822 & rec, gtype, vsize, &
823 & lbi, ubi, lbj, ubj, 1, n(ng), &
824 & scale, fmin, fmax, &
825# ifdef MASKING
826 & vmask, &
827# endif
828 & s_v(:,:,:,lout))
829 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
830 IF (master) THEN
831 WRITE (stdout,30) trim(vname(1,idvvel)), rec, trim(ncname)
832 END IF
833 exit_flag=3
834 ioerror=status
835 RETURN
836 ELSE
837 IF (master.and.lreport) THEN
838 WRITE (stdout,40) trim(vname(2,idvvel)), fmin, fmax
839 END IF
840 END IF
841
842# ifdef ADJUST_BOUNDARY
843!
844! Read in 3D V-momentum component open boundaries.
845!
846 IF (any(lobc(:,isvvel,ng))) THEN
847 ifield=idsbry(isvvel)
848 gtype=v3dvar
849 scale=1.0_dp
850 CALL netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
851 & my_ncid, varid)
852 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
853
854 status=nf_fread3d_bry(ng, model, ncname, my_ncid, &
855 & vname(1,ifield), varid, &
856 & rec, gtype, &
857 & lbij, ubij, 1, n(ng), nbrec(ng), &
858 & scale, fmin, fmax, &
859 & s_v_obc(:,:,:,:,lout))
860 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
861 IF (master) THEN
862 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
863 END IF
864 exit_flag=3
865 ioerror=status
866 RETURN
867 ELSE
868 IF (master.and.lreport) THEN
869 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
870 END IF
871 END IF
872 END IF
873# endif
874!
875! Read in tracers.
876!
877 gtype=r3dvar
878 scale=1.0_dp
879 DO itrc=1,nt(ng)
880 CALL netcdf_inq_varid (ng, model, ncname, &
881 & vname(1,idtvar(itrc)), my_ncid, varid)
882 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
883
884 status=nf_fread3d(ng, model, ncname, my_ncid, &
885 & vname(1,idtvar(itrc)), varid, &
886 & rec, gtype, vsize, &
887 & lbi, ubi, lbj, ubj, 1, n(ng), &
888 & scale, fmin, fmax, &
889# ifdef MASKING
890 & rmask, &
891# endif
892 & s_t(:,:,:,lout,itrc))
893 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
894 IF (master) THEN
895 WRITE (stdout,30) trim(vname(1,idtvar(itrc))), rec, &
896 & trim(ncname)
897 END IF
898 exit_flag=3
899 ioerror=status
900 RETURN
901 ELSE
902 IF (master.and.lreport) THEN
903 WRITE (stdout,40) trim(vname(2,idtvar(itrc))), fmin, fmax
904 END IF
905 END IF
906 END DO
907
908# ifdef ADJUST_BOUNDARY
909!
910! Read in tracers open boundaries.
911!
912 DO itrc=1,nt(ng)
913 IF (any(lobc(:,istvar(itrc),ng))) THEN
914 ifield=idsbry(istvar(itrc))
915 gtype=r3dvar
916 scale=1.0_dp
917 CALL netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
918 & my_ncid, varid)
919 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
920
921 status=nf_fread3d_bry(ng, model, ncname, my_ncid, &
922 & vname(1,ifield), varid, &
923 & rec, gtype, &
924 & lbij, ubij, 1, n(ng), nbrec(ng), &
925 & scale, fmin, fmax, &
926 & s_t_obc(:,:,:,:,lout,itrc))
927 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
928 IF (master) THEN
929 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
930 END IF
931 exit_flag=3
932 ioerror=status
933 RETURN
934 ELSE
935 IF (master.and.lreport) THEN
936 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
937 END IF
938 END IF
939 END IF
940 END DO
941# endif
942
943# ifdef ADJUST_STFLUX
944!
945! Read in surface tracers flux.
946!
947 gtype=r3dvar
948 scale=1.0_dp
949 DO itrc=1,nt(ng)
950 IF (lstflux(itrc,ng)) THEN
951 CALL netcdf_inq_varid (ng, model, ncname, &
952 & vname(1,idtsur(itrc)), my_ncid, varid)
953 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
954
955 status=nf_fread3d(ng, model, ncname, my_ncid, &
956 & vname(1,idtsur(itrc)), varid, &
957 & rec, gtype, vsize, &
958 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
959 & scale, fmin, fmax, &
960# ifdef MASKING
961 & rmask, &
962# endif
963 & s_tflux(:,:,:,lout,itrc))
964 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
965 IF (master) THEN
966 WRITE (stdout,30) trim(vname(1,idtsur(itrc))), rec, &
967 & trim(ncname)
968 END IF
969 exit_flag=3
970 ioerror=status
971 RETURN
972 ELSE
973 IF (master.and.lreport) THEN
974 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
975 END IF
976 END IF
977 END IF
978 END DO
979# endif
980#endif
981!
982! Close current file.
983!
984 IF (nopen.gt.0) THEN
985 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
986 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
987 END IF
988
989#ifdef PROFILE
990!
991! Turn off time wall clock.
992!
993 CALL wclock_off (ng, model, 80, __line__, myfile)
994#endif
995!
996 10 FORMAT ('STATE_READ_NF90 - ',a,t75,a, &
997 & /,19x,'(Grid ',i2.2,', t = ',a,', File: ',a, &
998 & ', Rec=',i4.4,', Index=',i1,')', &
999 & /,19x,'Called from ''',a,'''')
1000 20 FORMAT (' STATE_READ_NF90 - unable to open NetCDF file: ',a)
1001 30 FORMAT (' STATE_READ_NF90 - error while reading variable: ',a,2x, &
1002 & 'at time record = ',i3,/,14x,'in NetCDF file: ',a)
1003 40 FORMAT (16x,'- ',a,/,19x,'(Min = ',1p,e15.8, &
1004 & ' Max = ',1p,e15.8,')')
1005!
1006 RETURN
subroutine, public netcdf_close(ng, model, ncid, ncname, lupdate)
subroutine, public netcdf_open(ng, model, ncname, omode, ncid)
subroutine, public netcdf_sync(ng, model, ncname, ncid)
subroutine, public netcdf_inq_varid(ng, model, ncname, myvarname, ncid, varid)
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_iounits::calledfrom, mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::idfsur, mod_ncparam::idsbry, mod_ncparam::idtime, mod_ncparam::idtsur, mod_ncparam::idtvar, mod_ncparam::idubar, mod_ncparam::idusms, mod_ncparam::iduvel, mod_ncparam::idvbar, mod_ncparam::idvsms, mod_ncparam::idvvel, mod_iounits::ioerror, mod_ncparam::isfsur, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::lobc, mod_scalars::lstflux, mod_parallel::master, mod_param::n, mod_scalars::nbrec, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_inq_varid(), mod_netcdf::netcdf_open(), mod_netcdf::netcdf_sync(), mod_scalars::nfrec, mod_scalars::noerror, mod_param::nt, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::rclock, mod_scalars::sec2day, mod_iounits::sourcefile, mod_iounits::stdout, dateclock_mod::time_string(), mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, mod_param::v3dvar, mod_ncparam::vname, wclock_off(), and wclock_on().

Referenced by state_read().

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

◆ state_read_pio()

subroutine, private state_read_mod::state_read_pio ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
logical, intent(in) lreport,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lout,
integer, intent(in) rec,
integer, intent(in) nopen,
type (file_desc_t), intent(inout) piofile,
character (len=*), intent(in) ncname,
ifdef,
masking )
private

Definition at line 1012 of file state_read.F.

1017 & rmask, umask, vmask, &
1018# endif
1019# ifdef ADJUST_BOUNDARY
1020# ifdef SOLVE3D
1021 & s_t_obc, s_u_obc, s_v_obc, &
1022# endif
1023 & s_ubar_obc, s_vbar_obc, &
1024 & s_zeta_obc, &
1025# endif
1026# ifdef ADJUST_WSTRESS
1027 & s_ustr, s_vstr, &
1028# endif
1029# ifdef SOLVE3D
1030# ifdef ADJUST_STFLUX
1031 & s_tflux, &
1032# endif
1033 & s_t, s_u, s_v, &
1034# else
1035 & s_ubar, s_vbar, &
1036# endif
1037 & s_zeta)
1038!***********************************************************************
1039!
1040 USE mod_pio_netcdf
1041!
1042! Imported variable declarations.
1043!
1044 logical, intent(in) :: Lreport
1045!
1046 integer, intent(in) :: ng, tile, model
1047 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
1048 integer, intent(in) :: Lout, rec, nopen
1049!
1050 TYPE (File_desc_t), intent(inout) :: pioFile
1051!
1052 character (len=*), intent(in) :: ncname
1053!
1054# ifdef ASSUMED_SHAPE
1055# ifdef MASKING
1056 real(r8), intent(in) :: rmask(LBi:,LBj:)
1057 real(r8), intent(in) :: umask(LBi:,LBj:)
1058 real(r8), intent(in) :: vmask(LBi:,LBj:)
1059# endif
1060# ifdef ADJUST_BOUNDARY
1061# ifdef SOLVE3D
1062 real(r8), intent(inout) :: s_t_obc(LBij:,:,:,:,:,:)
1063 real(r8), intent(inout) :: s_u_obc(LBij:,:,:,:,:)
1064 real(r8), intent(inout) :: s_v_obc(LBij:,:,:,:,:)
1065# endif
1066 real(r8), intent(inout) :: s_ubar_obc(LBij:,:,:,:)
1067 real(r8), intent(inout) :: s_vbar_obc(LBij:,:,:,:)
1068 real(r8), intent(inout) :: s_zeta_obc(LBij:,:,:,:)
1069# endif
1070# ifdef ADJUST_WSTRESS
1071 real(r8), intent(inout) :: s_ustr(LBi:,LBj:,:,:)
1072 real(r8), intent(inout) :: s_vstr(LBi:,LBj:,:,:)
1073# endif
1074# ifdef SOLVE3D
1075# ifdef ADJUST_STFLUX
1076 real(r8), intent(inout) :: s_tflux(LBi:,LBj:,:,:,:)
1077# endif
1078 real(r8), intent(inout) :: s_t(LBi:,LBj:,:,:,:)
1079 real(r8), intent(inout) :: s_u(LBi:,LBj:,:,:)
1080 real(r8), intent(inout) :: s_v(LBi:,LBj:,:,:)
1081# else
1082 real(r8), intent(inout) :: s_ubar(LBi:,LBj:,:)
1083 real(r8), intent(inout) :: s_vbar(LBi:,LBj:,:)
1084# endif
1085 real(r8), intent(inout) :: s_zeta(LBi:,LBj:,:)
1086
1087# else
1088
1089# ifdef MASKING
1090 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
1091 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
1092 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
1093# endif
1094# ifdef ADJUST_BOUNDARY
1095# ifdef SOLVE3D
1096 real(r8), intent(inout) :: s_t_obc(LBij:UBij,N(ng),4, &
1097 & Nbrec(ng),2,NT(ng))
1098 real(r8), intent(inout) :: s_u_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
1099 real(r8), intent(inout) :: s_v_obc(LBij:UBij,N(ng),4,Nbrec(ng),2)
1100# endif
1101 real(r8), intent(inout) :: s_ubar_obc(LBij:UBij,4,Nbrec(ng),2)
1102 real(r8), intent(inout) :: s_vbar_obc(LBij:UBij,4,Nbrec(ng),2)
1103 real(r8), intent(inout) :: s_zeta_obc(LBij:UBij,4,Nbrec(ng),2)
1104# endif
1105# ifdef ADJUST_WSTRESS
1106 real(r8), intent(inout) :: s_ustr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
1107 real(r8), intent(inout) :: s_vstr(LBi:UBi,LBj:UBj,Nfrec(ng),2)
1108# endif
1109# ifdef SOLVE3D
1110# ifdef ADJUST_STFLUX
1111 real(r8), intent(inout) :: s_tflux(LBi:UBi,LBj:UBj, &
1112 & Nfrec(ng),2,NT(ng))
1113# endif
1114 real(r8), intent(inout) :: s_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
1115 real(r8), intent(inout) :: s_u(LBi:UBi,LBj:UBj,N(ng),2)
1116 real(r8), intent(inout) :: s_v(LBi:UBi,LBj:UBj,N(ng),2)
1117# else
1118 real(r8), intent(inout) :: s_ubar(LBi:UBi,LBj:UBj,:)
1119 real(r8), intent(inout) :: s_vbar(LBi:UBi,LBj:UBj,:)
1120# endif
1121 real(r8), intent(inout) :: s_zeta(LBi:UBi,LBj:UBj,:)
1122# endif
1123!
1124! Local variable declarations.
1125!
1126 integer :: Sstr, Send
1127 integer :: i, j, k
1128 integer :: ifield, itrc
1129 integer :: status
1130
1131 integer, dimension(4) :: Vsize
1132!
1133 real(r8) :: Fmin, Fmax
1134 real(dp) :: stime, scale
1135!
1136 character (len=15) :: Tstring
1137 character (len=22) :: t_code
1138
1139 character (len=*), parameter :: MyFile = &
1140 & __FILE__//", state_read_pio"
1141!
1142 TYPE (IO_desc_t), pointer :: my_ioDesc
1143 TYPE (File_desc_t) :: my_pioFile
1144 TYPE (My_VarDesc) :: my_pioVar
1145
1146# include "set_bounds.h"
1147!
1148 sourcefile=myfile
1149!
1150!-----------------------------------------------------------------------
1151! Read in requested model state record. Load data into state array
1152! index Lout.
1153!-----------------------------------------------------------------------
1154# ifdef PROFILE
1155!
1156! Turn on time wall clock.
1157!
1158 CALL wclock_on (ng, model, 80, __line__, myfile)
1159# endif
1160!
1161! Determine file and variables ids.
1162!
1163 IF ((nopen.gt.0).or.(piofile%fh.eq.-1)) THEN
1164 CALL pio_netcdf_open (ng, model, ncname, 0, my_piofile)
1165 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1166 piofile=my_piofile
1167 ELSE
1168 my_piofile=piofile
1169 END IF
1170!
1171 DO i=1,4
1172 vsize(i)=0
1173 END DO
1174
1175#ifdef SP4DVAR
1176!
1177! Synchronize NetCDF data to disk. It allows processes open for reading
1178! data to update the number of records written by other processes.
1179! The writer flushes buffers to disk, and the reader makes sure that
1180! it is reading from disk rather than from previously cached buffers.
1181!
1182 CALL pio_netcdf_sync (ng, model, ncname, my_piofile)
1183 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1184#endif
1185!
1186! Read in time.
1187!
1188 CALL pio_netcdf_get_time (ng, model, ncname, vname(1,idtime), &
1189 & rclock%DateNumber, stime, &
1190 & my_piofile, (/rec/), (/1/))
1191 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1192!
1193! Report information.
1194!
1195 IF (master) THEN
1196 CALL time_string (stime, t_code)
1197 sstr=scan(calledfrom,'/',back=.true.)+1
1198 send=len_trim(calledfrom)
1199 WRITE (tstring,'(f15.4)') stime*sec2day
1200 WRITE (stdout,10) 'Reading state fields,', t_code, &
1201 & ng, tstring, trim(ncname), rec, lout, &
1202 & calledfrom(sstr:send)
1203 END IF
1204!
1205! Read in free-surface.
1206!
1207 scale=1.0_dp
1208 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,idfsur), &
1209 & my_piofile, my_piovar%vd)
1210 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1211!
1212 my_piovar%gtype=r2dvar
1213 IF (kind(s_zeta).eq.8) THEN
1214 my_piovar%dkind=pio_double
1215 my_iodesc => iodesc_dp_r2dvar(ng)
1216 ELSE
1217 my_piovar%dkind=pio_real
1218 my_iodesc => iodesc_sp_r2dvar(ng)
1219 END IF
1220!
1221 status=nf_fread2d(ng, model, ncname, my_piofile, &
1222 & vname(1,idfsur), my_piovar, &
1223 & rec, my_iodesc, vsize, &
1224 & lbi, ubi, lbj, ubj, &
1225 & scale, fmin, fmax, &
1226# ifdef MASKING
1227 & rmask, &
1228# endif
1229 & s_zeta(:,:,lout))
1230 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1231 IF (master) THEN
1232 WRITE (stdout,30) trim(vname(1,idfsur)), rec, trim(ncname)
1233 END IF
1234 exit_flag=3
1235 ioerror=status
1236 RETURN
1237 ELSE
1238 IF (master.and.lreport) THEN
1239 WRITE (stdout,40) trim(vname(2,idfsur)), fmin, fmax
1240 END IF
1241 END IF
1242
1243# ifdef ADJUST_BOUNDARY
1244!
1245! Read in free-surface open boundaries.
1246!
1247 IF (any(lobc(:,isfsur,ng))) THEN
1248 ifield=idsbry(isfsur)
1249 scale=1.0_dp
1250 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
1251 & my_piofile, my_piovar%vd)
1252 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1253!
1254 my_piovar%gtype=r2dobc
1255 IF (kind(s_zeta_obc).eq.8) THEN
1256 my_piovar%dkind=pio_double
1257 my_iodesc => iodesc_dp_r2dobc(ng)
1258 ELSE
1259 my_piovar%dkind=pio_real
1260 my_iodesc => iodesc_sp_r2dobc(ng)
1261 END IF
1262!
1263 status=nf_fread2d_bry(ng, model, ncname, my_piofile, &
1264 & vname(1,ifield), my_piovar, &
1265 & rec, my_iodesc, &
1266 & lbij, ubij, nbrec(ng), &
1267 & scale, fmin, fmax, &
1268 & s_zeta_obc(:,:,:,lout))
1269 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1270 IF (master) THEN
1271 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
1272 END IF
1273 exit_flag=3
1274 ioerror=status
1275 RETURN
1276 ELSE
1277 IF (master.and.lreport) THEN
1278 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
1279 END IF
1280 END IF
1281 END IF
1282# endif
1283
1284# ifndef SOLVE3D
1285!
1286! Read in 2D U-momentum component.
1287!
1288 scale=1.0_dp
1289 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,idubar), &
1290 & my_piofile, my_piovar%vd)
1291 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1292!
1293 my_piovar%gtype=u2dvar
1294 IF (kind(s_ubar).eq.8) THEN
1295 my_piovar%dkind=pio_double
1296 my_iodesc => iodesc_dp_u2dvar(ng)
1297 ELSE
1298 my_piovar%dkind=pio_real
1299 my_iodesc => iodesc_sp_u2dvar(ng)
1300 END IF
1301!
1302 status=nf_fread2d(ng, model, ncname, my_piofile, &
1303 & vname(1,idubar), my_piovar, &
1304 & rec, my_iodesc, vsize, &
1305 & lbi, ubi, lbj, ubj, &
1306 & scale, fmin, fmax, &
1307# ifdef MASKING
1308 & umask, &
1309# endif
1310 & s_ubar(:,:,lout))
1311 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1312 IF (master) THEN
1313 WRITE (stdout,30) trim(vname(1,idubar)), rec, trim(ncname)
1314 END IF
1315 exit_flag=3
1316 ioerror=status
1317 RETURN
1318 ELSE
1319 IF (master.and.lreport) THEN
1320 WRITE (stdout,40) trim(vname(2,idubar)), fmin, fmax
1321 END IF
1322 END IF
1323!
1324! Read in 2D V-momentum component.
1325!
1326 scale=1.0_dp
1327 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,idvbar), &
1328 & my_piofile, my_piovar%vd)
1329 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1330!
1331 my_piovar%gtype=v2dvar
1332 IF (kind(s_vbar).eq.8) THEN
1333 my_piovar%dkind=pio_double
1334 my_iodesc => iodesc_dp_v2dvar(ng)
1335 ELSE
1336 my_piovar%dkind=pio_real
1337 my_iodesc => iodesc_sp_v2dvar(ng)
1338 END IF
1339!
1340 status=nf_fread2d(ng, model, ncname, my_piofile, &
1341 & vname(1,idvbar), my_piovar, &
1342 & rec, my_iodesc, vsize, &
1343 & lbi, ubi, lbj, ubj, &
1344 & scale, fmin, fmax, &
1345# ifdef MASKING
1346 & vmask, &
1347# endif
1348 & s_vbar(:,:,lout))
1349 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1350 IF (master) THEN
1351 WRITE (stdout,30) trim(vname(1,idvbar)), rec, trim(ncname)
1352 END IF
1353 exit_flag=3
1354 ioerror=status
1355 RETURN
1356 ELSE
1357 IF (master.and.lreport) THEN
1358 WRITE (stdout,40) trim(vname(2,idvbar)), fmin, fmax
1359 END IF
1360 END IF
1361# endif
1362
1363# ifdef ADJUST_BOUNDARY
1364!
1365! Read in 2D U-momentum component open boundaries.
1366!
1367 IF (any(lobc(:,isubar,ng))) THEN
1368 ifield=idsbry(isubar)
1369 scale=1.0_dp
1370 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
1371 & my_piofile, my_piovar%vd)
1372 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1373!
1374 my_piovar%gtype=u2dobc
1375 IF (kind(s_ubar_obc).eq.8) THEN
1376 my_piovar%dkind=pio_double
1377 my_iodesc => iodesc_dp_u2dobc(ng)
1378 ELSE
1379 my_piovar%dkind=pio_real
1380 my_iodesc => iodesc_sp_u2dobc(ng)
1381 END IF
1382!
1383 status=nf_fread2d_bry(ng, model, ncname, my_piofile, &
1384 & vname(1,ifield), my_piovar, &
1385 & rec, my_iodesc, &
1386 & lbij, ubij, nbrec(ng), &
1387 & scale, fmin, fmax, &
1388 & s_ubar_obc(:,:,:,lout))
1389 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1390 IF (master) THEN
1391 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
1392 END IF
1393 exit_flag=3
1394 ioerror=status
1395 RETURN
1396 ELSE
1397 IF (master.and.lreport) THEN
1398 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
1399 END IF
1400 END IF
1401 END IF
1402!
1403! Read in 2D V-momentum component open boundaries.
1404!
1405 IF (any(lobc(:,isvbar,ng))) THEN
1406 ifield=idsbry(isvbar)
1407 scale=1.0_dp
1408 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
1409 & my_piofile, my_piovar%vd)
1410 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1411!
1412 my_piovar%gtype=v2dobc
1413 IF (kind(s_vbar_obc).eq.8) THEN
1414 my_piovar%dkind=pio_double
1415 my_iodesc => iodesc_dp_v2dobc(ng)
1416 ELSE
1417 my_piovar%dkind=pio_real
1418 my_iodesc => iodesc_sp_v2dobc(ng)
1419 END IF
1420!
1421 status=nf_fread2d_bry(ng, model, ncname, my_piofile, &
1422 & vname(1,ifield), my_piovar, &
1423 & rec, my_iodesc, &
1424 & lbij, ubij, nbrec(ng), &
1425 & scale, fmin, fmax, &
1426 & s_vbar_obc(:,:,:,lout))
1427 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1428 IF (master) THEN
1429 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
1430 END IF
1431 exit_flag=3
1432 ioerror=status
1433 RETURN
1434 ELSE
1435 IF (master.and.lreport) THEN
1436 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
1437 END IF
1438 END IF
1439 END IF
1440# endif
1441
1442# ifdef ADJUST_WSTRESS
1443!
1444! Read surface momentum stress.
1445!
1446 scale=1.0_dp
1447 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,idusms), &
1448 & my_piofile, my_piovar%vd)
1449 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1450!
1451 my_piovar%gtype=u2dvar
1452 IF (kind(s_ustr).eq.8) THEN
1453 my_piovar%dkind=pio_double
1454 my_iodesc => iodesc_dp_u2dfrc(ng)
1455 ELSE
1456 my_piovar%dkind=pio_real
1457 my_iodesc => iodesc_sp_u2dfrc(ng)
1458 END IF
1459!
1460 status=nf_fread3d(ng, model, ncname, my_piofile, &
1461 & vname(1,idusms), my_piovar, &
1462 & rec, my_iodesc, vsize, &
1463 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
1464 & scale, fmin, fmax, &
1465# ifdef MASKING
1466 & umask, &
1467# endif
1468 & s_ustr(:,:,:,lout))
1469 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1470 IF (master) THEN
1471 WRITE (stdout,30) trim(vname(1,idusms)), rec, trim(ncname)
1472 END IF
1473 exit_flag=3
1474 ioerror=status
1475 RETURN
1476 ELSE
1477 IF (master.and.lreport) THEN
1478 WRITE (stdout,40) trim(vname(2,idusms)), fmin, fmax
1479 END IF
1480 END IF
1481!
1482 scale=1.0_dp
1483 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,idvsms), &
1484 & my_piofile, my_piovar%vd)
1485 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1486!
1487 my_piovar%gtype=v2dvar
1488 IF (kind(s_vstr).eq.8) THEN
1489 my_piovar%dkind=pio_double
1490 my_iodesc => iodesc_dp_v2dfrc(ng)
1491 ELSE
1492 my_piovar%dkind=pio_real
1493 my_iodesc => iodesc_sp_v2dfrc(ng)
1494 END IF
1495!
1496 status=nf_fread3d(ng, model, ncname, my_piofile, &
1497 & vname(1,idvsms), my_piovar, &
1498 & rec, my_iodesc, vsize, &
1499 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
1500 & scale, fmin, fmax, &
1501# ifdef MASKING
1502 & vmask, &
1503# endif
1504 & s_vstr(:,:,:,lout))
1505 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1506 IF (master) THEN
1507 WRITE (stdout,30) trim(vname(1,idvsms)), rec, trim(ncname)
1508 END IF
1509 exit_flag=3
1510 ioerror=status
1511 RETURN
1512 ELSE
1513 IF (master.and.lreport) THEN
1514 WRITE (stdout,40) trim(vname(2,idvsms)), fmin, fmax
1515 END IF
1516 END IF
1517# endif
1518
1519# ifdef SOLVE3D
1520!
1521! Read in 3D U-momentum component.
1522!
1523 scale=1.0_dp
1524 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,iduvel), &
1525 & my_piofile, my_piovar%vd)
1526 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1527!
1528 my_piovar%gtype=u3dvar
1529 IF (kind(s_u).eq.8) THEN
1530 my_piovar%dkind=pio_double
1531 my_iodesc => iodesc_dp_u3dvar(ng)
1532 ELSE
1533 my_piovar%dkind=pio_real
1534 my_iodesc => iodesc_sp_u3dvar(ng)
1535 END IF
1536!
1537 status=nf_fread3d(ng, model, ncname, my_piofile, &
1538 & vname(1,iduvel), my_piovar, &
1539 & rec, my_iodesc, vsize, &
1540 & lbi, ubi, lbj, ubj, 1, n(ng), &
1541 & scale, fmin, fmax, &
1542# ifdef MASKING
1543 & umask, &
1544# endif
1545 & s_u(:,:,:,lout))
1546 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1547 IF (master) THEN
1548 WRITE (stdout,30) trim(vname(1,iduvel)), rec, trim(ncname)
1549 END IF
1550 exit_flag=3
1551 ioerror=status
1552 RETURN
1553 ELSE
1554 IF (master.and.lreport) THEN
1555 WRITE (stdout,40) trim(vname(2,iduvel)), fmin, fmax
1556 END IF
1557 END IF
1558
1559# ifdef ADJUST_BOUNDARY
1560!
1561! Read in 3D U-momentum component open boundaries.
1562!
1563 IF (any(lobc(:,isuvel,ng))) THEN
1564 ifield=idsbry(isuvel)
1565 scale=1.0_dp
1566 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
1567 & my_piofile, my_piovar%vd)
1568 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1569!
1570 my_piovar%gtype=u3dobc
1571 IF (kind(s_u_obc).eq.8) THEN
1572 my_piovar%dkind=pio_double
1573 my_iodesc => iodesc_dp_u3dobc(ng)
1574 ELSE
1575 my_piovar%dkind=pio_real
1576 my_iodesc => iodesc_sp_u3dobc(ng)
1577 END IF
1578!
1579 status=nf_fread3d_bry(ng, model, ncname, my_piofile, &
1580 & vname(1,ifield), my_piovar, &
1581 & rec, my_iodesc, &
1582 & lbij, ubij, 1, n(ng), nbrec(ng), &
1583 & scale, fmin, fmax, &
1584 & s_u_obc(:,:,:,:,lout))
1585 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1586 IF (master) THEN
1587 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
1588 END IF
1589 exit_flag=3
1590 ioerror=status
1591 RETURN
1592 ELSE
1593 IF (master.and.lreport) THEN
1594 WRITE (stdout,40) trim(vname(1,ifield)), fmin, fmax
1595 END IF
1596 END IF
1597 END IF
1598# endif
1599!
1600! Read in 3D V-momentum component.
1601!
1602 scale=1.0_dp
1603 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,idvvel), &
1604 & my_piofile, my_piovar%vd)
1605 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1606!
1607 my_piovar%gtype=v3dvar
1608 IF (kind(s_v).eq.8) THEN
1609 my_piovar%dkind=pio_double
1610 my_iodesc => iodesc_dp_v3dvar(ng)
1611 ELSE
1612 my_piovar%dkind=pio_real
1613 my_iodesc => iodesc_sp_v3dvar(ng)
1614 END IF
1615!
1616 status=nf_fread3d(ng, model, ncname, my_piofile, &
1617 & vname(1,idvvel), my_piovar, &
1618 & rec, my_iodesc, vsize, &
1619 & lbi, ubi, lbj, ubj, 1, n(ng), &
1620 & scale, fmin, fmax, &
1621# ifdef MASKING
1622 & vmask, &
1623# endif
1624 & s_v(:,:,:,lout))
1625 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1626 IF (master) THEN
1627 WRITE (stdout,30) trim(vname(1,idvvel)), rec, trim(ncname)
1628 END IF
1629 exit_flag=3
1630 ioerror=status
1631 RETURN
1632 ELSE
1633 IF (master.and.lreport) THEN
1634 WRITE (stdout,40) trim(vname(2,idvvel)), fmin, fmax
1635 END IF
1636 END IF
1637
1638# ifdef ADJUST_BOUNDARY
1639!
1640! Read in 3D V-momentum component open boundaries.
1641!
1642 IF (any(lobc(:,isvvel,ng))) THEN
1643 ifield=idsbry(isvvel)
1644 scale=1.0_dp
1645 CALL pio_netcdf_inq_varid (ng, model, ncname, vname(1,ifield), &
1646 & my_piofile, my_piovar%vd)
1647 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1648!
1649 my_piovar%gtype=v3dobc
1650 IF (kind(s_v_obc).eq.8) THEN
1651 my_piovar%dkind=pio_double
1652 my_iodesc => iodesc_dp_v3dobc(ng)
1653 ELSE
1654 my_piovar%dkind=pio_real
1655 my_iodesc => iodesc_sp_v3dobc(ng)
1656 END IF
1657!
1658 status=nf_fread3d_bry(ng, model, ncname, my_piofile, &
1659 & vname(1,ifield), my_piovar, &
1660 & rec, my_iodesc, &
1661 & lbij, ubij, 1, n(ng), nbrec(ng), &
1662 & scale, fmin, fmax, &
1663 & s_v_obc(:,:,:,:,lout))
1664 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1665 IF (master) THEN
1666 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
1667 END IF
1668 exit_flag=3
1669 ioerror=status
1670 RETURN
1671 ELSE
1672 IF (master.and.lreport) THEN
1673 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
1674 END IF
1675 END IF
1676 END IF
1677# endif
1678!
1679! Read in tracers.
1680!
1681 scale=1.0_dp
1682 DO itrc=1,nt(ng)
1683 CALL pio_netcdf_inq_varid (ng, model, ncname, &
1684 & vname(1,idtvar(itrc)), &
1685 & my_piofile, my_piovar%vd)
1686 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1687!
1688 my_piovar%gtype=r3dvar
1689 IF (kind(s_t).eq.8) THEN
1690 my_piovar%dkind=pio_double
1691 my_iodesc => iodesc_dp_r3dvar(ng)
1692 ELSE
1693 my_piovar%dkind=pio_real
1694 my_iodesc => iodesc_sp_r3dvar(ng)
1695 END IF
1696!
1697 status=nf_fread3d(ng, model, ncname, my_piofile, &
1698 & vname(1,idtvar(itrc)), my_piovar, &
1699 & rec, my_iodesc, vsize, &
1700 & lbi, ubi, lbj, ubj, 1, n(ng), &
1701 & scale, fmin, fmax, &
1702# ifdef MASKING
1703 & rmask, &
1704# endif
1705 & s_t(:,:,:,lout,itrc))
1706 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1707 IF (master) THEN
1708 WRITE (stdout,30) trim(vname(1,idtvar(itrc))), rec, &
1709 & trim(ncname)
1710 END IF
1711 exit_flag=3
1712 ioerror=status
1713 RETURN
1714 ELSE
1715 IF (master.and.lreport) THEN
1716 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
1717 END IF
1718 END IF
1719 END DO
1720
1721# ifdef ADJUST_BOUNDARY
1722!
1723! Read in tracers open boundaries.
1724!
1725 DO itrc=1,nt(ng)
1726 IF (any(lobc(:,istvar(itrc),ng))) THEN
1727 ifield=idsbry(istvar(itrc))
1728 scale=1.0_dp
1729 CALL pio_netcdf_inq_varid (ng, model, ncname, &
1730 & vname(1,ifield), &
1731 & my_piofile, my_piovar%vd)
1732 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1733!
1734 my_piovar%gtype=r3dobc
1735 IF (kind(s_t_obc).eq.8) THEN
1736 my_piovar%dkind=pio_double
1737 my_iodesc => iodesc_dp_r3dobc(ng)
1738 ELSE
1739 my_piovar%dkind=pio_real
1740 my_iodesc => iodesc_sp_r3dobc(ng)
1741 END IF
1742!
1743 status=nf_fread3d_bry(ng, model, ncname, my_piofile, &
1744 & vname(1,ifield), my_piovar, &
1745 & rec, my_iodesc, &
1746 & lbij, ubij, 1, n(ng), nbrec(ng), &
1747 & scale, fmin, fmax, &
1748 & s_t_obc(:,:,:,:,lout,itrc))
1749 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1750 IF (master) THEN
1751 WRITE (stdout,30) trim(vname(1,ifield)), rec, trim(ncname)
1752 END IF
1753 exit_flag=3
1754 ioerror=status
1755 RETURN
1756 ELSE
1757 IF (master.and.lreport) THEN
1758 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
1759 END IF
1760 END IF
1761 END IF
1762 END DO
1763# endif
1764
1765# ifdef ADJUST_STFLUX
1766!
1767! Read in surface tracers flux.
1768!
1769 scale=1.0_dp
1770 DO itrc=1,nt(ng)
1771 IF (lstflux(itrc,ng)) THEN
1772 CALL pio_netcdf_inq_varid (ng, model, ncname, &
1773 & vname(1,idtsur(itrc)), &
1774 & my_piofile, my_piovar%vd)
1775 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1776!
1777 my_piovar%gtype=r2dvar
1778 IF (kind(s_tflux).eq.8) THEN
1779 my_piovar%dkind=pio_double
1780 my_iodesc => iodesc_dp_r2dfrc(ng)
1781 ELSE
1782 my_piovar%dkind=pio_real
1783 my_iodesc => iodesc_sp_r2dfrc(ng)
1784 END IF
1785!
1786 status=nf_fread3d(ng, model, ncname, my_piofile, &
1787 & vname(1,idtsur(itrc)), my_piovar, &
1788 & rec, my_iodesc, vsize, &
1789 & lbi, ubi, lbj, ubj, 1, nfrec(ng), &
1790 & scale, fmin, fmax, &
1791# ifdef MASKING
1792 & rmask, &
1793# endif
1794 & s_tflux(:,:,:,lout,itrc))
1795 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1796 IF (master) THEN
1797 WRITE (stdout,30) trim(vname(1,idtsur(itrc))), rec, &
1798 & trim(ncname)
1799 END IF
1800 exit_flag=3
1801 ioerror=status
1802 RETURN
1803 ELSE
1804 IF (master.and.lreport) THEN
1805 WRITE (stdout,40) trim(vname(2,ifield)), fmin, fmax
1806 END IF
1807 END IF
1808 END IF
1809 END DO
1810# endif
1811# endif
1812!
1813! Close current file.
1814!
1815 IF (nopen.gt.0) THEN
1816 CALL pio_netcdf_close (ng, model, my_piofile, ncname, .false.)
1817 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1818 END IF
1819
1820# ifdef PROFILE
1821!
1822! Turn off time wall clock.
1823!
1824 CALL wclock_off (ng, model, 80, __line__, myfile)
1825# endif
1826!
1827 10 FORMAT ('STATE_READ_PIO - ',a,t75,a, &
1828 & /,19x,'(Grid ',i2.2,', t = ',a,', File: ',a, &
1829 & ', Rec=',i4.4,', Index=',i1,')', &
1830 & /,19x,'Called from ''',a,'''')
1831 20 FORMAT (' STATE_READ_PIO - unable to open NetCDF file: ',a)
1832 30 FORMAT (' STATE_READ_PIO - error while reading variable: ',a,2x, &
1833 & 'at time record = ',i3,/,14x,'in NetCDF file: ',a)
1834 40 FORMAT (16x,'- ',a,/,19x,'(Min = ',1p,e15.8, &
1835 & ' Max = ',1p,e15.8,')')
1836!
1837 RETURN
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dobc
type(io_desc_t), dimension(:), pointer iodesc_dp_r3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_u3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_v3dobc
subroutine, public pio_netcdf_sync(ng, model, ncname, piofile)
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_u3dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_r3dobc
type(io_desc_t), dimension(:), pointer iodesc_dp_u3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_u3dobc
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dobc
subroutine, public pio_netcdf_close(ng, model, piofile, ncname, lupdate)
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dobc
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_v2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_v3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dfrc
subroutine, public pio_netcdf_inq_varid(ng, model, ncname, myvarname, piofile, piovar)
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_v3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_r3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r3dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_v3dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dfrc

References mod_iounits::calledfrom, mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::idfsur, mod_ncparam::idsbry, mod_ncparam::idtime, mod_ncparam::idtsur, mod_ncparam::idtvar, mod_ncparam::idubar, mod_ncparam::idusms, mod_ncparam::iduvel, mod_ncparam::idvbar, mod_ncparam::idvsms, mod_ncparam::idvvel, mod_pio_netcdf::iodesc_dp_r2dfrc, mod_pio_netcdf::iodesc_dp_r2dobc, mod_pio_netcdf::iodesc_dp_r2dvar, mod_pio_netcdf::iodesc_dp_r3dobc, mod_pio_netcdf::iodesc_dp_r3dvar, mod_pio_netcdf::iodesc_dp_u2dfrc, mod_pio_netcdf::iodesc_dp_u2dobc, mod_pio_netcdf::iodesc_dp_u2dvar, mod_pio_netcdf::iodesc_dp_u3dobc, mod_pio_netcdf::iodesc_dp_u3dvar, mod_pio_netcdf::iodesc_dp_v2dfrc, mod_pio_netcdf::iodesc_dp_v2dobc, mod_pio_netcdf::iodesc_dp_v2dvar, mod_pio_netcdf::iodesc_dp_v3dobc, mod_pio_netcdf::iodesc_dp_v3dvar, mod_pio_netcdf::iodesc_sp_r2dfrc, mod_pio_netcdf::iodesc_sp_r2dobc, mod_pio_netcdf::iodesc_sp_r2dvar, mod_pio_netcdf::iodesc_sp_r3dobc, mod_pio_netcdf::iodesc_sp_r3dvar, mod_pio_netcdf::iodesc_sp_u2dfrc, mod_pio_netcdf::iodesc_sp_u2dobc, mod_pio_netcdf::iodesc_sp_u2dvar, mod_pio_netcdf::iodesc_sp_u3dobc, mod_pio_netcdf::iodesc_sp_u3dvar, mod_pio_netcdf::iodesc_sp_v2dfrc, mod_pio_netcdf::iodesc_sp_v2dobc, mod_pio_netcdf::iodesc_sp_v2dvar, mod_pio_netcdf::iodesc_sp_v3dobc, mod_pio_netcdf::iodesc_sp_v3dvar, mod_iounits::ioerror, mod_ncparam::isfsur, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::lobc, mod_scalars::lstflux, mod_parallel::master, mod_param::n, mod_scalars::nbrec, mod_scalars::nfrec, mod_scalars::noerror, mod_param::nt, mod_pio_netcdf::pio_netcdf_close(), mod_pio_netcdf::pio_netcdf_inq_varid(), mod_pio_netcdf::pio_netcdf_open(), mod_pio_netcdf::pio_netcdf_sync(), mod_param::r2dobc, mod_param::r2dvar, mod_param::r3dobc, mod_param::r3dvar, mod_scalars::rclock, mod_scalars::sec2day, mod_iounits::sourcefile, mod_iounits::stdout, dateclock_mod::time_string(), mod_param::u2dobc, mod_param::u2dvar, mod_param::u3dobc, mod_param::u3dvar, mod_param::v2dobc, mod_param::v2dvar, mod_param::v3dobc, mod_param::v3dvar, mod_ncparam::vname, wclock_off(), and wclock_on().

Referenced by state_read().

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