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

Functions/Subroutines

subroutine, public time_corr (ng, tile, model, iotype, inpncname)
 
subroutine, private time_corr_nf90 (ng, tile, model, inpncname, lbi, ubi, lbj, ubj)
 
subroutine, private time_corr_pio (ng, tile, model, inpncname, lbi, ubi, lbj, ubj)
 

Function/Subroutine Documentation

◆ time_corr()

subroutine, public time_corr_mod::time_corr ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) iotype,
character (len=*), intent(in) inpncname )

Definition at line 57 of file time_corr.F.

58!***********************************************************************
59!
60! Imported variable declarations.
61!
62 integer, intent(in) :: ng, tile, model, IOtype
63!
64 character (len=*), intent(in) :: INPncname
65!
66! Local variable declarations.
67!
68 integer :: LBi, UBi, LBj, UBj
69!
70 character (len=*), parameter :: MyFile = &
71 & __FILE__
72!
73!-----------------------------------------------------------------------
74! Open existing nonlinear initial conditions NetCDF file and inquire
75! about its contains and define new variables, if needed.
76!-----------------------------------------------------------------------
77!
78 lbi=bounds(ng)%LBi(tile)
79 ubi=bounds(ng)%UBi(tile)
80 lbj=bounds(ng)%LBj(tile)
81 ubj=bounds(ng)%UBj(tile)
82!
83 SELECT CASE (iotype)
84 CASE (io_nf90)
85 CALL time_corr_nf90 (ng, tile, model, inpncname, &
86 & lbi, ubi, lbj, ubj)
87
88# if defined PIO_LIB && defined DISTRIBUTE
89 CASE (io_pio)
90 CALL time_corr_pio (ng, tile, model, inpncname, &
91 & lbi, ubi, lbj, ubj)
92# endif
93 CASE DEFAULT
94 IF (master) WRITE (stdout,10) iotype
95 exit_flag=2
96 END SELECT
97 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
98!
99 10 FORMAT (' TIME_CORR - Illegal input file type, io_type = ',i0, &
100 & /,13x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
101!
102 RETURN

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

Referenced by convolve_mod::error_covariance().

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

◆ time_corr_nf90()

subroutine, private time_corr_mod::time_corr_nf90 ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
character (len=*), intent(in) inpncname,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 106 of file time_corr.F.

108!***********************************************************************
109!
110 USE mod_netcdf
111!
112! Imported variable declarations.
113!
114 integer, intent(in) :: ng, tile, model
115 integer, intent(in) :: LBi, UBi, LBj, UBj
116!
117 character (len=*), intent(in) :: INPncname
118!
119! Local variable declarations.
120!
121 integer :: Iinp, Iout, Irec, MyRec, MyType, Nrec
122 integer :: INPncid, INPvid
123 integer :: i, gtype, status, varid
124 integer :: Vsize(4)
125 integer :: j, k, itrc
126!
127 real(r8) :: Fmin, Fmax
128!
129 real(dp) :: scale
130 real(dp) :: inp_time(1)
131 real(dp) :: timeI, timeR, fac
132!
133 character (len=*), parameter :: MyFile = &
134 & __FILE__//", time_corr_nf90"
135
136# include "set_bounds.h"
137!
138 sourcefile=myfile
139!
140!-----------------------------------------------------------------------
141! Determine variables to read and their availability.
142!-----------------------------------------------------------------------
143!
144! Inquire about the dimensions and check for consistency.
145!
146 CALL netcdf_check_dim (ng, model, inpncname)
147 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
148 nrec=rec_size
149!
150! Inquire about the variables.
151!
152 CALL netcdf_inq_var (ng, model, inpncname)
153 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
154!
155! Set Vsize to zero to deactivate interpolation of input data to model
156! grid in "nf_fread2d".
157!
158 DO i=1,4
159 vsize(i)=0
160 END DO
161!
162!-----------------------------------------------------------------------
163! Read adjoint solution for time convolutions.
164!-----------------------------------------------------------------------
165!
166! Open input NetCDF file.
167!
168 CALL netcdf_open (ng, model, inpncname, 0, inpncid)
169 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
170 WRITE (stdout,10) trim(inpncname)
171 RETURN
172 END IF
173!
174! Compute the time convolutions for the model error corrections.
175!
176 iinp=1
177 scale=1.0_dp
178!
179! Read all free-surface records into the "f_zetaS" array.
180!
181 IF (find_string(var_name, n_var, vname(1,idztlf), inpvid)) THEN
182 gtype=var_flag(inpvid)*r2dvar
183 mytype=gtype
184 DO myrec=1,nrec-1
185 status=nf_fread2d(ng, model, inpncname, inpncid, &
186 & vname(1,idztlf), inpvid, &
187 & myrec, mytype, vsize, &
188 & lbi, ubi, lbj, ubj, &
189 & scale, fmin, fmax, &
190# ifdef MASKING
191 & grid(ng) % rmask, &
192# endif
193 & ocean(ng) % f_zetaS(:,:,myrec))
194 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
195 IF (master) THEN
196 WRITE (stdout,20) trim(vname(1,idztlf)), myrec, &
197 & trim(inpncname)
198 END IF
199 exit_flag=2
200 ioerror=status
201 RETURN
202 END IF
203 END DO
204 ELSE
205 IF (master) WRITE (stdout,40) trim(vname(1,idztlf)), &
206 & trim(inpncname)
207 exit_flag=2
208 RETURN
209 END IF
210!
211! Read all 2D U-momentum records into the "f_ubarS" array.
212!
213 IF (find_string(var_name, n_var, vname(1,idubtf), inpvid)) THEN
214 gtype=var_flag(inpvid)*u2dvar
215 mytype=gtype
216 DO myrec=1,nrec-1
217 status=nf_fread2d(ng, model, inpncname, inpncid, &
218 & vname(1,idubtf), inpvid, &
219 & myrec, mytype, vsize, &
220 & lbi, ubi, lbj, ubj, &
221 & scale, fmin, fmax, &
222# ifdef MASKING
223 & grid(ng) % umask_full, &
224# endif
225 & ocean(ng) % f_ubarS(:,:,myrec))
226 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
227 IF (master) THEN
228 WRITE (stdout,20) trim(vname(1,idubtf)), myrec, &
229 & trim(inpncname)
230 END IF
231 exit_flag=2
232 ioerror=status
233 RETURN
234 END IF
235 END DO
236 ELSE
237 IF (master) WRITE (stdout,20) trim(vname(1,idubtf)), &
238 & trim(inpncname)
239 exit_flag=2
240 RETURN
241 END IF
242!
243! Read all 2D V-momentum records into the "f_vbarS" array.
244!
245 IF (find_string(var_name, n_var, vname(1,idvbtf), inpvid)) THEN
246 gtype=var_flag(inpvid)*v2dvar
247 mytype=gtype
248 DO myrec=1,nrec-1
249 status=nf_fread2d(ng, model, inpncname, inpncid, &
250 & vname(1,idvbtf), inpvid, &
251 & myrec, mytype, vsize, &
252 & lbi, ubi, lbj, ubj, &
253 & scale, fmin, fmax, &
254# ifdef MASKING
255 & grid(ng) % vmask_full, &
256# endif
257 & ocean(ng) % f_vbarS(:,:,myrec))
258 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
259 IF (master) THEN
260 WRITE (stdout,20) trim(vname(1,idvbtf)), myrec, &
261 & trim(inpncname)
262 END IF
263 exit_flag=2
264 ioerror=status
265 RETURN
266 END IF
267 END DO
268 ELSE
269 IF (master) WRITE (stdout,40) trim(vname(1,idvbtf)), &
270 & trim(inpncname)
271 exit_flag=2
272 RETURN
273 END IF
274
275# ifdef SOLVE3D
276!
277! Read all 3D U-momentum records into the "f_uS" array.
278!
279 IF (find_string(var_name, n_var, vname(1,idutlf), inpvid)) THEN
280 gtype=var_flag(inpvid)*u3dvar
281 mytype=gtype
282 DO myrec=1,nrec-1
283 status=nf_fread3d(ng, model, inpncname, inpncid, &
284 & vname(1,idutlf), inpvid, &
285 & myrec, mytype, vsize, &
286 & lbi, ubi, lbj, ubj, 1, n(ng), &
287 & scale, fmin, fmax, &
288# ifdef MASKING
289 & grid(ng) % umask_full, &
290# endif
291 & ocean(ng) % f_uS(:,:,:,myrec))
292 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
293 IF (master) THEN
294 WRITE (stdout,20) trim(vname(1,idutlf)), myrec, &
295 & trim(inpncname)
296 END IF
297 exit_flag=2
298 ioerror=status
299 RETURN
300 END IF
301 END DO
302 ELSE
303 IF (master) WRITE (stdout,40) trim(vname(1,idutlf)), &
304 & trim(inpncname)
305 exit_flag=2
306 RETURN
307 END IF
308!
309! Read all 3D V-momentum records into the "f_vS" array.
310!
311 IF (find_string(var_name, n_var, vname(1,idvtlf), inpvid)) THEN
312 gtype=var_flag(inpvid)*v3dvar
313 mytype=gtype
314 DO myrec=1,nrec-1
315 status=nf_fread3d(ng, model, inpncname, inpncid, &
316 & vname(1,idvtlf), inpvid, &
317 & myrec, mytype, vsize, &
318 & lbi, ubi, lbj, ubj, 1, n(ng), &
319 & scale, fmin, fmax, &
320# ifdef MASKING
321 & grid(ng) % vmask_full, &
322# endif
323 & ocean(ng) % f_vS(:,:,:,myrec))
324 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
325 IF (master) THEN
326 WRITE (stdout,20) trim(vname(1,idvtlf)), myrec, &
327 & trim(inpncname)
328 END IF
329 exit_flag=2
330 ioerror=status
331 RETURN
332 END IF
333 END DO
334 ELSE
335 IF (master) WRITE (stdout,40) trim(vname(1,idvtlf)), &
336 & trim(inpncname)
337 exit_flag=2
338 RETURN
339 END IF
340!
341! Read all the tracers records into the "f_tS" array.
342!
343 DO itrc=1,nt(ng)
344 IF (find_string(var_name, n_var, vname(1,idttlf(itrc)), &
345 & inpvid)) THEN
346 gtype=var_flag(inpvid)*r3dvar
347 mytype=gtype
348 DO myrec=1,nrec-1
349 status=nf_fread3d(ng, model, inpncname, inpncid, &
350 & vname(1,idttlf(itrc)), inpvid, &
351 & myrec, mytype, vsize, &
352 & lbi, ubi, lbj, ubj, 1, n(ng), &
353 & scale, fmin, fmax, &
354# ifdef MASKING
355 & grid(ng) % rmask, &
356# endif
357 & ocean(ng) % f_tS(:,:,:,myrec,itrc))
358 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
359 IF (master) THEN
360 WRITE (stdout,20) trim(vname(1,idttlf(itrc))), myrec, &
361 & trim(inpncname)
362 END IF
363 exit_flag=2
364 ioerror=status
365 RETURN
366 END IF
367 END DO
368 ELSE
369 IF (master) WRITE (stdout,40) trim(vname(1,idttlf(i))), &
370 & trim(inpncname)
371 exit_flag=2
372 RETURN
373 END IF
374 END DO
375# endif
376!
377!=======================================================================
378! Perform time convolution and write impulse forcing into output TLF
379! NetCDF file.
380!=======================================================================
381!
382 iout=0
383 rec_loop : DO irec=nrec-1,1,-1
384 iout=iout+1
385 timei=dstart*day2sec+(iout-1)*nadj(ng)*dt(ng)
386!
387! Write out time variable.
388!
389 IF (find_string(var_name, n_var, vname(1,idtime), inpvid)) THEN
390 CALL netcdf_get_time (ng, model, inpncname, vname(1,idtime), &
391 & rclock%DateNumber, inp_time, &
392 & ncid = inpncid, &
393 & start = (/irec/), total = (/1/))
394 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
395
396 CALL netcdf_put_fvar (ng, model, tlf(ng)%name, &
397 & vname(1,idtime), inp_time, &
398 & (/iout/), (/1/), &
399 & ncid = tlf(ng)%ncid)
400 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
401 ELSE
402 IF (master) WRITE (stdout,20) trim(vname(1,idtime)), &
403 & trim(inpncname)
404 exit_flag=2
405 END IF
406!
407!-----------------------------------------------------------------------
408! Free-surface weak-constraint time convolution.
409!-----------------------------------------------------------------------
410!
411! Clear adjoint state array.
412!
413 DO j=jstrr,jendr
414 DO i=istrr,iendr
415 ocean(ng) % ad_zeta(i,j,iinp)=0.0_r8
416 END DO
417 END DO
418!
419! Convolve.
420!
421 IF (find_string(var_name, n_var, vname(1,idztlf), inpvid)) THEN
422 gtype=var_flag(inpvid)*r2dvar
423 mytype=gtype
424 DO myrec=nrec-1,1,-1
425!
426! Compute weighting factors.
427!
428 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
429 fac=exp(-abs(timei-timer)/tdecay(isfsur,ng))
430# ifdef ENDPOINT_TRAPEZOIDAL
431 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
432 fac=0.5_r8*fac
433 END IF
434# endif
435!
436! Form the weighted sum of the adjoint file records in time.
437!
438 DO j=jstrr,jendr
439 DO i=istrr,iendr
440 ocean(ng) % ad_zeta(i,j,iinp)= &
441 & ocean(ng) % ad_zeta(i,j,iinp)+ &
442 & fac*ocean(ng) % f_zetaS(i,j,myrec)
443 END DO
444 END DO
445 END DO
446!
447! Write out convolved solution.
448!
449 mytype=gtype
450 status=nf_fwrite2d(ng, model, tlf(ng)%ncid, idztlf, &
451 & tlf(ng)%Vid(idztlf), &
452 & iout, mytype, &
453 & lbi, ubi, lbj, ubj, scale, &
454# ifdef MASKING
455 & grid(ng) % rmask, &
456# endif
457 & ocean(ng) % ad_zeta(:,:,iinp))
458 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
459 IF (master) THEN
460 WRITE (stdout,30) trim(vname(1,idztlf)), irec, &
461 & trim(tlf(ng)%name)
462 END IF
463 exit_flag=3
464 ioerror=status
465 RETURN
466 END IF
467 ELSE
468 IF (master) WRITE (stdout,40) trim(vname(1,idztlf)), &
469 & trim(inpncname)
470 exit_flag=2
471 RETURN
472 END IF
473!
474!-----------------------------------------------------------------------
475! 2D U-momentum weak-constraint time convolution.
476!-----------------------------------------------------------------------
477!
478! Clear adjoint state array.
479!
480 DO j=jstrr,jendr
481 DO i=istr,iendr
482 ocean(ng) % ad_ubar(i,j,iinp)=0.0_r8
483 END DO
484 END DO
485!
486! Convolve.
487!
488 IF (find_string(var_name, n_var, vname(1,idubtf), inpvid)) THEN
489 gtype=var_flag(inpvid)*u2dvar
490 mytype=gtype
491 DO myrec=nrec-1,1,-1
492!
493! Compute weighting factors.
494!
495 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
496 fac=exp(-abs(timei-timer)/tdecay(isubar,ng))
497# ifdef ENDPOINT_TRAPEZOIDAL
498 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
499 fac=0.5_r8*fac
500 END IF
501# endif
502!
503! Form the weighted sum of the adjoint file records in time.
504!
505 DO j=jstrr,jendr
506 DO i=istr,iendr
507 ocean(ng) % ad_ubar(i,j,iinp)= &
508 & ocean(ng) % ad_ubar(i,j,iinp)+ &
509 & fac*ocean(ng) % f_ubarS(i,j,myrec)
510 END DO
511 END DO
512 END DO
513!
514! Write out convolved solution.
515!
516 mytype=gtype
517 status=nf_fwrite2d(ng, model, tlf(ng)%ncid, idubtf, &
518 & tlf(ng)%Vid(idubtf), &
519 & iout, mytype, &
520 & lbi, ubi, lbj, ubj, scale, &
521# ifdef MASKING
522 & grid(ng) % umask_full, &
523# endif
524 & ocean(ng) % ad_ubar(:,:,iinp))
525 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
526 IF (master) THEN
527 WRITE (stdout,30) trim(vname(1,idubtf)), irec, &
528 & trim(tlf(ng)%name)
529 END IF
530 exit_flag=3
531 ioerror=status
532 RETURN
533 END IF
534 ELSE
535 IF (master) WRITE (stdout,20) trim(vname(1,idubtf)), &
536 & trim(inpncname)
537 exit_flag=2
538 RETURN
539 END IF
540!
541!-----------------------------------------------------------------------
542! 2D V-momentum weak-constraint time convolution.
543!-----------------------------------------------------------------------
544!
545! Clear adjoint state array.
546!
547 DO j=jstr,jendr
548 DO i=istrr,iendr
549 ocean(ng) % ad_vbar(i,j,iinp)=0.0_r8
550 END DO
551 END DO
552!
553! Convolve.
554!
555 IF (find_string(var_name, n_var, vname(1,idvbtf), inpvid)) THEN
556 gtype=var_flag(inpvid)*v2dvar
557 mytype=gtype
558 DO myrec=nrec-1,1,-1
559!
560! Compute weighting factors.
561!
562 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
563 fac=exp(-abs(timei-timer)/tdecay(isvbar,ng))
564# ifdef ENDPOINT_TRAPEZOIDAL
565 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
566 fac=0.5_r8*fac
567 END IF
568# endif
569!
570! Form the weighted sum of the adjoint file records in time.
571!
572 DO j=jstr,jendr
573 DO i=istrr,iendr
574 ocean(ng) % ad_vbar(i,j,iinp)= &
575 & ocean(ng) % ad_vbar(i,j,iinp)+ &
576 & fac*ocean(ng) % f_vbarS(i,j,myrec)
577 END DO
578 END DO
579 END DO
580!
581! Write out convolved solution.
582!
583 status=nf_fwrite2d(ng, model, tlf(ng)%ncid, idvbtf, &
584 & tlf(ng)%Vid(idvbtf), &
585 & iout, mytype, &
586 & lbi, ubi, lbj, ubj, scale, &
587# ifdef MASKING
588 & grid(ng) % vmask_full, &
589# endif
590 & ocean(ng) % ad_vbar(:,:,iinp))
591 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
592 IF (master) THEN
593 WRITE (stdout,30) trim(vname(1,idvbtf)), irec, &
594 & trim(tlf(ng)%name)
595 END IF
596 exit_flag=3
597 ioerror=status
598 RETURN
599 END IF
600 ELSE
601 IF (master) WRITE (stdout,20) trim(vname(1,idvbtf)), &
602 & trim(inpncname)
603 exit_flag=2
604 RETURN
605 END IF
606
607# ifdef SOLVE3D
608!
609!-----------------------------------------------------------------------
610! 3D U-momentum weak-constraint time convolution.
611!-----------------------------------------------------------------------
612!
613! Clear adjoint state array.
614!
615 DO k=1,n(ng)
616 DO j=jstrr,jendr
617 DO i=istr,iendr
618 ocean(ng) % ad_u(i,j,k,iinp)=0.0_r8
619 END DO
620 END DO
621 END DO
622!
623! Convolve.
624!
625 IF (find_string(var_name, n_var, vname(1,idutlf), inpvid)) THEN
626 gtype=var_flag(inpvid)*u3dvar
627 mytype=gtype
628 DO myrec=nrec-1,1,-1
629!
630! Compute weighting factors.
631!
632 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
633 fac=exp(-abs(timei-timer)/tdecay(isuvel,ng))
634# ifdef ENDPOINT_TRAPEZOIDAL
635 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
636 fac=0.5_r8*fac
637 END IF
638# endif
639!
640! Form the weighted sum of the adjoint file records in time.
641!
642 DO k=1,n(ng)
643 DO j=jstrr,jendr
644 DO i=istr,iendr
645 ocean(ng) % ad_u(i,j,k,iinp)= &
646 & ocean(ng) % ad_u(i,j,k,iinp)+ &
647 & fac*ocean(ng) % f_uS(i,j,k,myrec)
648 END DO
649 END DO
650 END DO
651 END DO
652!
653! Write out convolved solution.
654!
655 status=nf_fwrite3d(ng, model, tlf(ng)%ncid, idutlf, &
656 & tlf(ng)%Vid(idutlf), &
657 & iout, mytype, &
658 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
659# ifdef MASKING
660 & grid(ng) % umask_full, &
661# endif
662 & ocean(ng) % ad_u(:,:,:,iinp))
663 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
664 IF (master) THEN
665 WRITE (stdout,30) trim(vname(1,idutlf)), irec, &
666 & trim(tlf(ng)%name)
667 END IF
668 exit_flag=3
669 ioerror=status
670 RETURN
671 END IF
672 ELSE
673 IF (master) WRITE (stdout,40) trim(vname(1,idutlf)), &
674 & trim(inpncname)
675 exit_flag=2
676 RETURN
677 END IF
678!
679!-----------------------------------------------------------------------
680! 3D V-momentum weak-constraint time convolution.
681!-----------------------------------------------------------------------
682!
683! Clear adjoint state array.
684!
685 DO k=1,n(ng)
686 DO j=jstr,jendr
687 DO i=istrr,iendr
688 ocean(ng) % ad_v(i,j,k,iinp)=0.0_r8
689 END DO
690 END DO
691 END DO
692!
693! Convolve.
694!
695 IF (find_string(var_name, n_var, vname(1,idvtlf), inpvid)) THEN
696 gtype=var_flag(inpvid)*v3dvar
697 mytype=gtype
698 DO myrec=nrec-1,1,-1
699!
700! Compute weighting factors.
701!
702 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
703 fac=exp(-abs(timei-timer)/tdecay(isvvel,ng))
704# ifdef ENDPOINT_TRAPEZOIDAL
705 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
706 fac=0.5_r8*fac
707 END IF
708# endif
709!
710! Form the weighted sum of the adjoint file records in time.
711!
712 DO k=1,n(ng)
713 DO j=jstr,jendr
714 DO i=istrr,iendr
715 ocean(ng) % ad_v(i,j,k,iinp)= &
716 & ocean(ng) % ad_v(i,j,k,iinp)+ &
717 & fac*ocean(ng) % f_vS(i,j,k,myrec)
718 END DO
719 END DO
720 END DO
721 END DO
722!
723! Write out convolved solution.
724!
725 status=nf_fwrite3d(ng, model, tlf(ng)%ncid, idvtlf, &
726 & tlf(ng)%Vid(idvtlf), &
727 & iout, mytype, &
728 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
729# ifdef MASKING
730 & grid(ng) % vmask_full, &
731# endif
732 & ocean(ng) % ad_v(:,:,:,iinp))
733 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
734 IF (master) THEN
735 WRITE (stdout,30) trim(vname(1,idvtlf)), irec, &
736 & trim(tlf(ng)%name)
737 END IF
738 exit_flag=3
739 ioerror=status
740 RETURN
741 END IF
742 ELSE
743 IF (master) WRITE (stdout,40) trim(vname(1,idvtlf)), &
744 & trim(inpncname)
745 exit_flag=2
746 RETURN
747 END IF
748!
749!-----------------------------------------------------------------------
750! Tracers type variables weak-constraint time convolution.
751!-----------------------------------------------------------------------
752!
753! Clear adjoint state array.
754!
755 DO itrc=1,nt(ng)
756 DO k=1,n(ng)
757 DO j=jstrr,jendr
758 DO i=istrr,iendr
759 ocean(ng) % ad_t(i,j,k,iinp,itrc)=0.0_r8
760 END DO
761 END DO
762 END DO
763 END DO
764!
765! Convolve.
766!
767 DO itrc=1,nt(ng)
768 IF (find_string(var_name, n_var, vname(1,idttlf(itrc)), &
769 & inpvid)) THEN
770 gtype=var_flag(inpvid)*r3dvar
771 mytype=gtype
772 DO myrec=nrec-1,1,-1
773!
774! Compute weighting factors.
775!
776 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
777 fac=exp(-abs(timei-timer)/tdecay(istvar(itrc),ng))
778# ifdef ENDPOINT_TRAPEZOIDAL
779 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
780 fac=0.5_r8*fac
781 END IF
782# endif
783!
784! Form the weighted sum of the adjoint file records in time.
785!
786 DO k=1,n(ng)
787 DO j=jstr,jendr
788 DO i=istrr,iendr
789 ocean(ng) % ad_t(i,j,k,iinp,itrc)= &
790 & ocean(ng) % ad_t(i,j,k,iinp,itrc)+ &
791 & fac*ocean(ng) % f_tS(i,j,k,myrec,itrc)
792 END DO
793 END DO
794 END DO
795 END DO
796!
797! Write out convolved solution.
798!
799 status=nf_fwrite3d(ng, model, tlf(ng)%ncid, idttlf(itrc), &
800 & tlf(ng)%Tid(itrc), &
801 & iout, mytype, &
802 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
803# ifdef MASKING
804 & grid(ng) % rmask, &
805# endif
806 & ocean(ng) % ad_t(:,:,:,iinp,itrc))
807 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
808 IF (master) THEN
809 WRITE (stdout,30) trim(vname(1,idttlf(itrc))), irec, &
810 & trim(tlf(ng)%name)
811 END IF
812 exit_flag=3
813 ioerror=status
814 RETURN
815 END IF
816 ELSE
817 IF (master) WRITE (stdout,40) trim(vname(1,idttlf(i))), &
818 & trim(inpncname)
819 exit_flag=2
820 RETURN
821 END IF
822 END DO
823# endif
824 END DO rec_loop
825!
826!-----------------------------------------------------------------------
827! Synchronize impulse NetCDF file to disk to allow other processes
828! to access data immediately after it is written.
829!-----------------------------------------------------------------------
830!
831 CALL netcdf_sync (ng, model, tlf(ng)%name, tlf(ng)%ncid)
832 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
833!
834 IF (master) WRITE (stdout,50) nrec-1, trim(tlf(ng)%name)
835!
836 10 FORMAT (/,' TIME_CORR_NF90 - unable to open input NetCDF', &
837 & ' file: ',a)
838 20 FORMAT (/,' TIME_CORR_NF90 - error while reading variable: ',a, &
839 & 2x,'at time record = ',i3,/,18x,'in input NetCDF file:', &
840 & 1x,a)
841 30 FORMAT (/,' TIME_CORR_NF90 - error while writing variable: ',a, &
842 & 2x,'at time record = ',i3,/,18x,'into NetCDF file: ',a)
843 40 FORMAT (/,' TIME_CORR_NF90 - cannot find state variable: ',a, &
844 & /,18x,'in input NetCDF file: ',a)
845 50 FORMAT (2x,'TIME_CORR_NF90 - wrote convolved adjoint impulses', &
846 & ', records: 001 to ',i3.3,/,21x,'file: ',a)
847!
848 RETURN
integer, dimension(mvars) var_flag
Definition mod_netcdf.F:162
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_sync(ng, model, ncname, ncid)
integer rec_size
Definition mod_netcdf.F:156
subroutine, public netcdf_inq_var(ng, model, ncname, ncid, myvarname, searchvar, varid, nvardim, nvaratt)

References mod_scalars::day2sec, mod_scalars::dstart, mod_scalars::dt, mod_scalars::exit_flag, strings_mod::find_string(), strings_mod::founderror(), mod_grid::grid, mod_ncparam::idtime, mod_ncparam::idttlf, mod_ncparam::idubtf, mod_ncparam::idutlf, mod_ncparam::idvbtf, mod_ncparam::idvtlf, mod_ncparam::idztlf, mod_iounits::ioerror, mod_ncparam::isfsur, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_parallel::master, mod_param::n, mod_netcdf::n_var, mod_scalars::nadj, mod_netcdf::netcdf_check_dim(), mod_netcdf::netcdf_inq_var(), mod_netcdf::netcdf_open(), mod_netcdf::netcdf_sync(), mod_scalars::noerror, mod_param::nt, mod_ocean::ocean, mod_param::r2dvar, mod_param::r3dvar, mod_scalars::rclock, mod_netcdf::rec_size, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::tdecay, mod_iounits::tlf, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, mod_param::v3dvar, mod_netcdf::var_flag, mod_netcdf::var_name, and mod_ncparam::vname.

Referenced by time_corr().

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

◆ time_corr_pio()

subroutine, private time_corr_mod::time_corr_pio ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
character (len=*), intent(in) inpncname,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 854 of file time_corr.F.

856!***********************************************************************
857!
859!
860! Imported variable declarations.
861!
862 integer, intent(in) :: ng, tile, model
863 integer, intent(in) :: LBi, UBi, LBj, UBj
864!
865 character (len=*), intent(in) :: INPncname
866!
867! Local variable declarations.
868!
869 integer :: Iinp, Iout, Irec, MyRec, Nrec
870 integer :: i, status, vindex
871 integer :: Vsize(4)
872 integer :: j, k, itrc
873!
874 real(r8) :: Fmin, Fmax
875!
876 real(dp) :: scale
877 real(dp) :: inp_time(1)
878 real(dp) :: timeI, timeR, fac
879!
880 character (len=*), parameter :: MyFile = &
881 & __FILE__//", time_corr_pio"
882!
883 TYPE (IO_Desc_t), pointer :: ioDesc
884 TYPE (My_VarDesc), pointer :: pioVar
885
886# include "set_bounds.h"
887!
888 sourcefile=myfile
889!
890!-----------------------------------------------------------------------
891! Determine variables to read and their availability.
892!-----------------------------------------------------------------------
893!
894! Inquire about the dimensions and check for consistency.
895!
896 CALL pio_netcdf_check_dim (ng, model, inpncname)
897 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
898 nrec=rec_size
899!
900! Inquire about the variables.
901!
902 CALL pio_netcdf_inq_var (ng, model, inpncname)
903 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
904!
905! Set Vsize to zero to deactivate interpolation of input data to model
906! grid in "nf_fread2d".
907!
908 DO i=1,4
909 vsize(i)=0
910 END DO
911!
912!-----------------------------------------------------------------------
913! Read adjoint solution for time convolutions.
914!-----------------------------------------------------------------------
915!
916! Open input NetCDF file.
917!
918 CALL pio_netcdf_open (ng, model, inpncname, 0, inppiofile)
919 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
920 WRITE (stdout,10) trim(inpncname)
921 RETURN
922 END IF
923!
924! Compute the time convolutions for the model error corrections.
925!
926 iinp=1
927 scale=1.0_dp
928!
929! Read all free-surface records into the "f_zetaS" array.
930!
931 IF (find_string(var_name, n_var, vname(1,idztlf), vindex)) THEN
932
933 piovar%vd => var_desc(vindex)
934 IF (kind(ocean(ng)%f_zetaS).eq.8) THEN
935 piovar%dkind=pio_double
936 iodesc => iodesc_dp_r2dvar(ng)
937 ELSE
938 piovar%dkind=pio_real
939 iodesc => iodesc_sp_r2dvar(ng)
940 END IF
941 piovar%gtype=r2dvar
942!
943 DO myrec=1,nrec-1
944 status=nf_fread2d(ng, model, inpncname, inppiofile, &
945 & vname(1,idztlf), piovar, &
946 & myrec, iodesc, vsize, &
947 & lbi, ubi, lbj, ubj, &
948 & scale, fmin, fmax, &
949# ifdef MASKING
950 & grid(ng) % rmask, &
951# endif
952 & ocean(ng) % f_zetaS(:,:,myrec))
953 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
954 IF (master) THEN
955 WRITE (stdout,20) trim(vname(1,idztlf)), myrec, &
956 & trim(inpncname)
957 END IF
958 exit_flag=2
959 ioerror=status
960 RETURN
961 END IF
962 END DO
963 ELSE
964 IF (master) WRITE (stdout,40) trim(vname(1,idztlf)), &
965 & trim(inpncname)
966 exit_flag=2
967 RETURN
968 END IF
969!
970! Read all 2D U-momentum records into the "f_ubarS" array.
971!
972 IF (find_string(var_name, n_var, vname(1,idubtf), vindex)) THEN
973
974 piovar%vd => var_desc(vindex)
975 IF (kind(ocean(ng)%f_ubarS).eq.8) THEN
976 piovar%dkind=pio_double
977 iodesc => iodesc_dp_u2dvar(ng)
978 ELSE
979 piovar%dkind=pio_real
980 iodesc => iodesc_sp_u2dvar(ng)
981 END IF
982 piovar%gtype=u2dvar
983!
984 DO myrec=1,nrec-1
985 status=nf_fread2d(ng, model, inpncname, inppiofile, &
986 & vname(1,idubtf), piovar, &
987 & myrec, iodesc, vsize, &
988 & lbi, ubi, lbj, ubj, &
989 & scale, fmin, fmax, &
990# ifdef MASKING
991 & grid(ng) % umask_full, &
992# endif
993 & ocean(ng) % f_ubarS(:,:,myrec))
994 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
995 IF (master) THEN
996 WRITE (stdout,20) trim(vname(1,idubtf)), myrec, &
997 & trim(inpncname)
998 END IF
999 exit_flag=2
1000 ioerror=status
1001 RETURN
1002 END IF
1003 END DO
1004 ELSE
1005 IF (master) WRITE (stdout,20) trim(vname(1,idubtf)), &
1006 & trim(inpncname)
1007 exit_flag=2
1008 RETURN
1009 END IF
1010!
1011! Read all 2D V-momentum records into the "f_vbarS" array.
1012!
1013 IF (find_string(var_name, n_var, vname(1,idvbtf), vindex)) THEN
1014
1015 piovar%vd => var_desc(vindex)
1016 IF (kind(ocean(ng)%f_vbarS).eq.8) THEN
1017 piovar%dkind=pio_double
1018 iodesc => iodesc_dp_v2dvar(ng)
1019 ELSE
1020 piovar%dkind=pio_real
1021 iodesc => iodesc_sp_v2dvar(ng)
1022 END IF
1023 piovar%gtype=v2dvar
1024!
1025 DO myrec=1,nrec-1
1026 status=nf_fread2d(ng, model, inpncname, inppiofile, &
1027 & vname(1,idvbtf), piovar, &
1028 & myrec, iodesc, vsize, &
1029 & lbi, ubi, lbj, ubj, &
1030 & scale, fmin, fmax, &
1031# ifdef MASKING
1032 & grid(ng) % vmask_full, &
1033# endif
1034 & ocean(ng) % f_vbarS(:,:,myrec))
1035 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1036 IF (master) THEN
1037 WRITE (stdout,20) trim(vname(1,idvbtf)), myrec, &
1038 & trim(inpncname)
1039 END IF
1040 exit_flag=2
1041 ioerror=status
1042 RETURN
1043 END IF
1044 END DO
1045 ELSE
1046 IF (master) WRITE (stdout,40) trim(vname(1,idvbtf)), &
1047 & trim(inpncname)
1048 exit_flag=2
1049 RETURN
1050 END IF
1051
1052# ifdef SOLVE3D
1053!
1054! Read all 3D U-momentum records into the "f_uS" array.
1055!
1056 IF (find_string(var_name, n_var, vname(1,idutlf), vindex)) THEN
1057
1058 piovar%vd => var_desc(vindex)
1059 IF (kind(ocean(ng)%f_uS).eq.8) THEN
1060 piovar%dkind=pio_double
1061 iodesc => iodesc_dp_u3dvar(ng)
1062 ELSE
1063 piovar%dkind=pio_real
1064 iodesc => iodesc_sp_u3dvar(ng)
1065 END IF
1066 piovar%gtype=u3dvar
1067!
1068 DO myrec=1,nrec-1
1069 status=nf_fread3d(ng, model, inpncname, inppiofile, &
1070 & vname(1,idutlf), piovar, &
1071 & myrec, iodesc, vsize, &
1072 & lbi, ubi, lbj, ubj, 1, n(ng), &
1073 & scale, fmin, fmax, &
1074# ifdef MASKING
1075 & grid(ng) % umask_full, &
1076# endif
1077 & ocean(ng) % f_uS(:,:,:,myrec))
1078 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1079 IF (master) THEN
1080 WRITE (stdout,20) trim(vname(1,idutlf)), myrec, &
1081 & trim(inpncname)
1082 END IF
1083 exit_flag=2
1084 ioerror=status
1085 RETURN
1086 END IF
1087 END DO
1088 ELSE
1089 IF (master) WRITE (stdout,40) trim(vname(1,idutlf)), &
1090 & trim(inpncname)
1091 exit_flag=2
1092 RETURN
1093 END IF
1094!
1095! Read all 3D V-momentum records into the "f_vS" array.
1096!
1097 IF (find_string(var_name, n_var, vname(1,idvtlf), vindex)) THEN
1098
1099 piovar%vd => var_desc(vindex)
1100 IF (kind(ocean(ng)%f_vS).eq.8) THEN
1101 piovar%dkind=pio_double
1102 iodesc => iodesc_dp_v3dvar(ng)
1103 ELSE
1104 piovar%dkind=pio_real
1105 iodesc => iodesc_sp_v3dvar(ng)
1106 END IF
1107 piovar%gtype=v3dvar
1108!
1109 DO myrec=1,nrec-1
1110 status=nf_fread3d(ng, model, inpncname, inppiofile, &
1111 & vname(1,idvtlf), piovar, &
1112 & myrec, iodesc, vsize, &
1113 & lbi, ubi, lbj, ubj, 1, n(ng), &
1114 & scale, fmin, fmax, &
1115# ifdef MASKING
1116 & grid(ng) % vmask_full, &
1117# endif
1118 & ocean(ng) % f_vS(:,:,:,myrec))
1119 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1120 IF (master) THEN
1121 WRITE (stdout,20) trim(vname(1,idvtlf)), myrec, &
1122 & trim(inpncname)
1123 END IF
1124 exit_flag=2
1125 ioerror=status
1126 RETURN
1127 END IF
1128 END DO
1129 ELSE
1130 IF (master) WRITE (stdout,40) trim(vname(1,idvtlf)), &
1131 & trim(inpncname)
1132 exit_flag=2
1133 RETURN
1134 END IF
1135!
1136! Read all the tracers records into the "f_tS" array.
1137!
1138 DO itrc=1,nt(ng)
1139 IF (find_string(var_name, n_var, vname(1,idttlf(itrc)), &
1140 & vindex)) THEN
1141
1142 piovar%vd => var_desc(vindex)
1143 IF (kind(ocean(ng)%f_tS).eq.8) THEN
1144 piovar%dkind=pio_double
1145 iodesc => iodesc_dp_r3dvar(ng)
1146 ELSE
1147 piovar%dkind=pio_real
1148 iodesc => iodesc_sp_r3dvar(ng)
1149 END IF
1150 piovar%gtype=r3dvar
1151!
1152 DO myrec=1,nrec-1
1153 status=nf_fread3d(ng, model, inpncname, inppiofile, &
1154 & vname(1,idttlf(itrc)), piovar, &
1155 & myrec, iodesc, vsize, &
1156 & lbi, ubi, lbj, ubj, 1, n(ng), &
1157 & scale, fmin, fmax, &
1158# ifdef MASKING
1159 & grid(ng) % rmask, &
1160# endif
1161 & ocean(ng) % f_tS(:,:,:,myrec,itrc))
1162 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1163 IF (master) THEN
1164 WRITE (stdout,20) trim(vname(1,idttlf(itrc))), myrec, &
1165 & trim(inpncname)
1166 END IF
1167 exit_flag=2
1168 ioerror=status
1169 RETURN
1170 END IF
1171 END DO
1172 ELSE
1173 IF (master) WRITE (stdout,40) trim(vname(1,idttlf(i))), &
1174 & trim(inpncname)
1175 exit_flag=2
1176 RETURN
1177 END IF
1178 END DO
1179# endif
1180!
1181!=======================================================================
1182! Perform time convolution and write impulse forcing into output TLF
1183! NetCDF file.
1184!=======================================================================
1185!
1186 iout=0
1187 rec_loop : DO irec=nrec-1,1,-1
1188 iout=iout+1
1189 timei=dstart*day2sec+(iout-1)*nadj(ng)*dt(ng)
1190!
1191! Write out time variable.
1192!
1193 IF (find_string(var_name, n_var, vname(1,idtime), vindex)) THEN
1194 CALL pio_netcdf_get_time (ng, model, inpncname, &
1195 & vname(1,idtime), &
1196 & rclock%DateNumber, inp_time, &
1197 & piofile = inppiofile, &
1198 & start = (/irec/), &
1199 & total = (/1/))
1200 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1201!
1202 CALL pio_netcdf_put_fvar (ng, model, tlf(ng)%name, &
1203 & vname(1,idtime), inp_time, &
1204 & (/iout/), (/1/), &
1205 & piofile = tlf(ng)%pioFile)
1206 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1207 ELSE
1208 IF (master) WRITE (stdout,20) trim(vname(1,idtime)), &
1209 & trim(inpncname)
1210 exit_flag=2
1211 END IF
1212!
1213!-----------------------------------------------------------------------
1214! Free-surface weak-constraint time convolution.
1215!-----------------------------------------------------------------------
1216!
1217! Clear adjoint state array.
1218!
1219 DO j=jstrr,jendr
1220 DO i=istrr,iendr
1221 ocean(ng) % ad_zeta(i,j,iinp)=0.0_r8
1222 END DO
1223 END DO
1224!
1225! Convolve.
1226!
1227 IF (find_string(var_name, n_var, vname(1,idztlf), vindex)) THEN
1228
1229 IF (tlf(ng)%pioVar(idztlf)%dkind.eq.pio_double) THEN
1230 iodesc => iodesc_dp_r2dvar(ng)
1231 ELSE
1232 iodesc => iodesc_sp_r2dvar(ng)
1233 END IF
1234
1235 DO myrec=nrec-1,1,-1
1236!
1237! Compute weighting factors.
1238!
1239 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
1240 fac=exp(-abs(timei-timer)/tdecay(isfsur,ng))
1241# ifdef ENDPOINT_TRAPEZOIDAL
1242 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
1243 fac=0.5_r8*fac
1244 END IF
1245# endif
1246!
1247! Form the weighted sum of the adjoint file records in time.
1248!
1249 DO j=jstrr,jendr
1250 DO i=istrr,iendr
1251 ocean(ng) % ad_zeta(i,j,iinp)= &
1252 & ocean(ng) % ad_zeta(i,j,iinp)+ &
1253 & fac*ocean(ng) % f_zetaS(i,j,myrec)
1254 END DO
1255 END DO
1256 END DO
1257!
1258! Write out convolved solution.
1259!
1260 status=nf_fwrite2d(ng, model, tlf(ng)%pioFile, idztlf, &
1261 & tlf(ng)%pioVar(idztlf), &
1262 & iout, iodesc, &
1263 & lbi, ubi, lbj, ubj, scale, &
1264# ifdef MASKING
1265 & grid(ng) % rmask, &
1266# endif
1267 & ocean(ng) % ad_zeta(:,:,iinp))
1268 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1269 IF (master) THEN
1270 WRITE (stdout,30) trim(vname(1,idztlf)), irec, &
1271 & trim(tlf(ng)%name)
1272 END IF
1273 exit_flag=3
1274 ioerror=status
1275 RETURN
1276 END IF
1277 ELSE
1278 IF (master) WRITE (stdout,40) trim(vname(1,idztlf)), &
1279 & trim(inpncname)
1280 exit_flag=2
1281 RETURN
1282 END IF
1283!
1284!-----------------------------------------------------------------------
1285! 2D U-momentum weak-constraint time convolution.
1286!-----------------------------------------------------------------------
1287!
1288! Clear adjoint state array.
1289!
1290 DO j=jstrr,jendr
1291 DO i=istr,iendr
1292 ocean(ng) % ad_ubar(i,j,iinp)=0.0_r8
1293 END DO
1294 END DO
1295!
1296! Convolve.
1297!
1298 IF (find_string(var_name, n_var, vname(1,idubtf), vindex)) THEN
1299
1300 IF (tlf(ng)%pioVar(idubtf)%dkind.eq.pio_double) THEN
1301 iodesc => iodesc_dp_u2dvar(ng)
1302 ELSE
1303 iodesc => iodesc_sp_u2dvar(ng)
1304 END IF
1305
1306 DO myrec=nrec-1,1,-1
1307!
1308! Compute weighting factors.
1309!
1310 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
1311 fac=exp(-abs(timei-timer)/tdecay(isubar,ng))
1312# ifdef ENDPOINT_TRAPEZOIDAL
1313 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
1314 fac=0.5_r8*fac
1315 END IF
1316# endif
1317!
1318! Form the weighted sum of the adjoint file records in time.
1319!
1320 DO j=jstrr,jendr
1321 DO i=istr,iendr
1322 ocean(ng) % ad_ubar(i,j,iinp)= &
1323 & ocean(ng) % ad_ubar(i,j,iinp)+ &
1324 & fac*ocean(ng) % f_ubarS(i,j,myrec)
1325 END DO
1326 END DO
1327 END DO
1328!
1329! Write out convolved solution.
1330!
1331 status=nf_fwrite2d(ng, model, tlf(ng)%pioFile, idubtf, &
1332 & tlf(ng)%pioVar(idubtf), &
1333 & iout, iodesc, &
1334 & lbi, ubi, lbj, ubj, scale, &
1335# ifdef MASKING
1336 & grid(ng) % umask_full, &
1337# endif
1338 & ocean(ng) % ad_ubar(:,:,iinp))
1339 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1340 IF (master) THEN
1341 WRITE (stdout,30) trim(vname(1,idubtf)), irec, &
1342 & trim(tlf(ng)%name)
1343 END IF
1344 exit_flag=3
1345 ioerror=status
1346 RETURN
1347 END IF
1348 ELSE
1349 IF (master) WRITE (stdout,20) trim(vname(1,idubtf)), &
1350 & trim(inpncname)
1351 exit_flag=2
1352 RETURN
1353 END IF
1354!
1355!-----------------------------------------------------------------------
1356! 2D V-momentum weak-constraint time convolution.
1357!-----------------------------------------------------------------------
1358!
1359! Clear adjoint state array.
1360!
1361 DO j=jstr,jendr
1362 DO i=istrr,iendr
1363 ocean(ng) % ad_vbar(i,j,iinp)=0.0_r8
1364 END DO
1365 END DO
1366!
1367! Convolve.
1368!
1369 IF (find_string(var_name, n_var, vname(1,idvbtf), vindex)) THEN
1370
1371 IF (tlf(ng)%pioVar(idvbtf)%dkind.eq.pio_double) THEN
1372 iodesc => iodesc_dp_v2dvar(ng)
1373 ELSE
1374 iodesc => iodesc_sp_v2dvar(ng)
1375 END IF
1376!
1377 DO myrec=nrec-1,1,-1
1378!
1379! Compute weighting factors.
1380!
1381 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
1382 fac=exp(-abs(timei-timer)/tdecay(isvbar,ng))
1383# ifdef ENDPOINT_TRAPEZOIDAL
1384 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
1385 fac=0.5_r8*fac
1386 END IF
1387# endif
1388!
1389! Form the weighted sum of the adjoint file records in time.
1390!
1391 DO j=jstr,jendr
1392 DO i=istrr,iendr
1393 ocean(ng) % ad_vbar(i,j,iinp)= &
1394 & ocean(ng) % ad_vbar(i,j,iinp)+ &
1395 & fac*ocean(ng) % f_vbarS(i,j,myrec)
1396 END DO
1397 END DO
1398 END DO
1399!
1400! Write out convolved solution.
1401!
1402 status=nf_fwrite2d(ng, model, tlf(ng)%pioFile, idvbtf, &
1403 & tlf(ng)%pioVar(idvbtf), &
1404 & iout, iodesc, &
1405 & lbi, ubi, lbj, ubj, scale, &
1406# ifdef MASKING
1407 & grid(ng) % vmask_full, &
1408# endif
1409 & ocean(ng) % ad_vbar(:,:,iinp))
1410 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1411 IF (master) THEN
1412 WRITE (stdout,30) trim(vname(1,idvbtf)), irec, &
1413 & trim(tlf(ng)%name)
1414 END IF
1415 exit_flag=3
1416 ioerror=status
1417 RETURN
1418 END IF
1419 ELSE
1420 IF (master) WRITE (stdout,20) trim(vname(1,idvbtf)), &
1421 & trim(inpncname)
1422 exit_flag=2
1423 RETURN
1424 END IF
1425
1426# ifdef SOLVE3D
1427!
1428!-----------------------------------------------------------------------
1429! 3D U-momentum weak-constraint time convolution.
1430!-----------------------------------------------------------------------
1431!
1432! Clear adjoint state array.
1433!
1434 DO k=1,n(ng)
1435 DO j=jstrr,jendr
1436 DO i=istr,iendr
1437 ocean(ng) % ad_u(i,j,k,iinp)=0.0_r8
1438 END DO
1439 END DO
1440 END DO
1441!
1442! Convolve.
1443!
1444 IF (find_string(var_name, n_var, vname(1,idutlf), vindex)) THEN
1445
1446 IF (tlf(ng)%pioVar(idutlf)%dkind.eq.pio_double) THEN
1447 iodesc => iodesc_dp_u3dvar(ng)
1448 ELSE
1449 iodesc => iodesc_sp_u3dvar(ng)
1450 END IF
1451!
1452 DO myrec=nrec-1,1,-1
1453!
1454! Compute weighting factors.
1455!
1456 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
1457 fac=exp(-abs(timei-timer)/tdecay(isuvel,ng))
1458# ifdef ENDPOINT_TRAPEZOIDAL
1459 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
1460 fac=0.5_r8*fac
1461 END IF
1462# endif
1463!
1464! Form the weighted sum of the adjoint file records in time.
1465!
1466 DO k=1,n(ng)
1467 DO j=jstrr,jendr
1468 DO i=istr,iendr
1469 ocean(ng) % ad_u(i,j,k,iinp)= &
1470 & ocean(ng) % ad_u(i,j,k,iinp)+ &
1471 & fac*ocean(ng) % f_uS(i,j,k,myrec)
1472 END DO
1473 END DO
1474 END DO
1475 END DO
1476!
1477! Write out convolved solution.
1478!
1479 status=nf_fwrite3d(ng, model, tlf(ng)%pioFile, idutlf, &
1480 & tlf(ng)%pioVar(idutlf), &
1481 & iout, iodesc, &
1482 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
1483# ifdef MASKING
1484 & grid(ng) % umask_full, &
1485# endif
1486 & ocean(ng) % ad_u(:,:,:,iinp))
1487 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1488 IF (master) THEN
1489 WRITE (stdout,30) trim(vname(1,idutlf)), irec, &
1490 & trim(tlf(ng)%name)
1491 END IF
1492 exit_flag=3
1493 ioerror=status
1494 RETURN
1495 END IF
1496 ELSE
1497 IF (master) WRITE (stdout,40) trim(vname(1,idutlf)), &
1498 & trim(inpncname)
1499 exit_flag=2
1500 RETURN
1501 END IF
1502!
1503!-----------------------------------------------------------------------
1504! 3D V-momentum weak-constraint time convolution.
1505!-----------------------------------------------------------------------
1506!
1507! Clear adjoint state array.
1508!
1509 DO k=1,n(ng)
1510 DO j=jstr,jendr
1511 DO i=istrr,iendr
1512 ocean(ng) % ad_v(i,j,k,iinp)=0.0_r8
1513 END DO
1514 END DO
1515 END DO
1516!
1517! Convolve.
1518!
1519 IF (find_string(var_name, n_var, vname(1,idvtlf), vindex)) THEN
1520
1521 IF (tlf(ng)%pioVar(idvtlf)%dkind.eq.pio_double) THEN
1522 iodesc => iodesc_dp_v3dvar(ng)
1523 ELSE
1524 iodesc => iodesc_sp_v3dvar(ng)
1525 END IF
1526!
1527 DO myrec=nrec-1,1,-1
1528!
1529! Compute weighting factors.
1530!
1531 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
1532 fac=exp(-abs(timei-timer)/tdecay(isvvel,ng))
1533# ifdef ENDPOINT_TRAPEZOIDAL
1534 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
1535 fac=0.5_r8*fac
1536 END IF
1537# endif
1538!
1539! Form the weighted sum of the adjoint file records in time.
1540!
1541 DO k=1,n(ng)
1542 DO j=jstr,jendr
1543 DO i=istrr,iendr
1544 ocean(ng) % ad_v(i,j,k,iinp)= &
1545 & ocean(ng) % ad_v(i,j,k,iinp)+ &
1546 & fac*ocean(ng) % f_vS(i,j,k,myrec)
1547 END DO
1548 END DO
1549 END DO
1550 END DO
1551!
1552! Write out convolved solution.
1553!
1554 status=nf_fwrite3d(ng, model, tlf(ng)%pioFile, idvtlf, &
1555 & tlf(ng)%pioVar(idvtlf), &
1556 & iout, iodesc, &
1557 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
1558# ifdef MASKING
1559 & grid(ng) % vmask_full, &
1560# endif
1561 & ocean(ng) % ad_v(:,:,:,iinp))
1562 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1563 IF (master) THEN
1564 WRITE (stdout,30) trim(vname(1,idvtlf)), irec, &
1565 & trim(tlf(ng)%name)
1566 END IF
1567 exit_flag=3
1568 ioerror=status
1569 RETURN
1570 END IF
1571 ELSE
1572 IF (master) WRITE (stdout,40) trim(vname(1,idvtlf)), &
1573 & trim(inpncname)
1574 exit_flag=2
1575 RETURN
1576 END IF
1577!
1578!-----------------------------------------------------------------------
1579! Tracers type variables weak-constraint time convolution.
1580!-----------------------------------------------------------------------
1581!
1582! Clear adjoint state array.
1583!
1584 DO itrc=1,nt(ng)
1585 DO k=1,n(ng)
1586 DO j=jstrr,jendr
1587 DO i=istrr,iendr
1588 ocean(ng) % ad_t(i,j,k,iinp,itrc)=0.0_r8
1589 END DO
1590 END DO
1591 END DO
1592 END DO
1593!
1594! Convolve.
1595!
1596 DO itrc=1,nt(ng)
1597 IF (find_string(var_name, n_var, vname(1,idttlf(itrc)), &
1598 & vindex)) THEN
1599
1600 IF (tlf(ng)%pioTrc(itrc)%dkind.eq.pio_double) THEN
1601 iodesc => iodesc_dp_r3dvar(ng)
1602 ELSE
1603 iodesc => iodesc_sp_r3dvar(ng)
1604 END IF
1605!
1606 DO myrec=nrec-1,1,-1
1607!
1608! Compute weighting factors.
1609!
1610 timer=dstart*day2sec+(nrec-1-myrec)*nadj(ng)*dt(ng)
1611 fac=exp(-abs(timei-timer)/tdecay(istvar(itrc),ng))
1612# ifdef ENDPOINT_TRAPEZOIDAL
1613 IF ((myrec.eq.1).or.(myrec.eq.(nrec-1))) THEN
1614 fac=0.5_r8*fac
1615 END IF
1616# endif
1617!
1618! Form the weighted sum of the adjoint file records in time.
1619!
1620 DO k=1,n(ng)
1621 DO j=jstr,jendr
1622 DO i=istrr,iendr
1623 ocean(ng) % ad_t(i,j,k,iinp,itrc)= &
1624 & ocean(ng) % ad_t(i,j,k,iinp,itrc)+ &
1625 & fac*ocean(ng) % f_tS(i,j,k,myrec,itrc)
1626 END DO
1627 END DO
1628 END DO
1629 END DO
1630!
1631! Write out convolved solution.
1632!
1633 status=nf_fwrite3d(ng, model, tlf(ng)%pioFile, &
1634 & idttlf(itrc), tlf(ng)%pioTrc(itrc), &
1635 & iout, iodesc, &
1636 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
1637# ifdef MASKING
1638 & grid(ng) % rmask, &
1639# endif
1640 & ocean(ng) % ad_t(:,:,:,iinp,itrc))
1641 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1642 IF (master) THEN
1643 WRITE (stdout,30) trim(vname(1,idttlf(itrc))), irec, &
1644 & trim(tlf(ng)%name)
1645 END IF
1646 exit_flag=3
1647 ioerror=status
1648 RETURN
1649 END IF
1650 ELSE
1651 IF (master) WRITE (stdout,40) trim(vname(1,idttlf(i))), &
1652 & trim(inpncname)
1653 exit_flag=2
1654 RETURN
1655 END IF
1656 END DO
1657# endif
1658 END DO rec_loop
1659!
1660!-----------------------------------------------------------------------
1661! Synchronize impulse NetCDF file to disk to allow other processes
1662! to access data immediately after it is written.
1663!-----------------------------------------------------------------------
1664!
1665 CALL pio_netcdf_sync (ng, model, tlf(ng)%name, tlf(ng)%pioFile)
1666 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1667
1668 IF (master) WRITE (stdout,50) nrec-1, trim(tlf(ng)%name)
1669!
1670 10 FORMAT (/,' TIME_CORR_PIO - unable to open input NetCDF', &
1671 & ' file: ',a)
1672 20 FORMAT (/,' TIME_CORR_PIO - error while reading variable: ',a, &
1673 & 2x,'at time record = ',i3,/,17x,'in input NetCDF file:', &
1674 & 1x,a)
1675 30 FORMAT (/,' TIME_CORR_PIO - error while writing variable: ',a, &
1676 & 2x,'at time record = ',i3,/,17x,'into NetCDF file: ',a)
1677 40 FORMAT (/,' TIME_CORR_PIO - cannot find state variable: ',a, &
1678 & /,18x,'in input NetCDF file: ',a)
1679 50 FORMAT (2x,'TIME_CORR_PIO - wrote convolved adjoint impulses', &
1680 & ', records: 001 to ',i3.3,/,21x,'file: ',a)
1681!
1682 RETURN
type(io_desc_t), dimension(:), pointer iodesc_dp_r3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_u3dvar
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_u3dvar
type(var_desc_t), dimension(:), pointer var_desc
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dvar
subroutine, public pio_netcdf_inq_var(ng, model, ncname, piofile, myvarname, searchvar, piovar, nvardim, nvaratt)
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_v3dvar
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_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_sp_v2dvar

References mod_scalars::day2sec, mod_scalars::dstart, mod_scalars::dt, mod_scalars::exit_flag, strings_mod::find_string(), strings_mod::founderror(), mod_grid::grid, mod_ncparam::idtime, mod_ncparam::idttlf, mod_ncparam::idubtf, mod_ncparam::idutlf, mod_ncparam::idvbtf, mod_ncparam::idvtlf, mod_ncparam::idztlf, mod_pio_netcdf::iodesc_dp_r2dvar, mod_pio_netcdf::iodesc_dp_r3dvar, mod_pio_netcdf::iodesc_dp_u2dvar, mod_pio_netcdf::iodesc_dp_u3dvar, mod_pio_netcdf::iodesc_dp_v2dvar, mod_pio_netcdf::iodesc_dp_v3dvar, mod_pio_netcdf::iodesc_sp_r2dvar, mod_pio_netcdf::iodesc_sp_r3dvar, mod_pio_netcdf::iodesc_sp_u2dvar, mod_pio_netcdf::iodesc_sp_u3dvar, mod_pio_netcdf::iodesc_sp_v2dvar, 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_parallel::master, mod_param::n, mod_scalars::nadj, mod_scalars::noerror, mod_param::nt, mod_ocean::ocean, mod_pio_netcdf::pio_netcdf_check_dim(), mod_pio_netcdf::pio_netcdf_inq_var(), mod_pio_netcdf::pio_netcdf_open(), mod_pio_netcdf::pio_netcdf_sync(), mod_param::r2dvar, mod_param::r3dvar, mod_scalars::rclock, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::tdecay, mod_iounits::tlf, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, mod_param::v3dvar, mod_pio_netcdf::var_desc, and mod_ncparam::vname.

Referenced by time_corr().

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