ROMS
Loading...
Searching...
No Matches
obs_read.F
Go to the documentation of this file.
1#include "cppdefs.h"
3#if (defined FOUR_DVAR || defined VERIFICATION) && defined OBSERVATIONS
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This module reads in observations data when appropriate from the !
13! input file using the standard NetCDF library or the Parallel-IO !
14! (PIO) library. !
15! !
16! The observations data is stored for use elsewhere. !
17! !
18!=======================================================================
19!
20 USE mod_param
21 USE mod_parallel
22 USE mod_fourdvar
23 USE mod_iounits
24 USE mod_ncparam
25 USE mod_scalars
26!
27 USE dateclock_mod, ONLY : time_string
28 USE strings_mod, ONLY : founderror
29!
30 implicit none
31!
32 PUBLIC :: obs_read
33 PRIVATE :: obs_read_nf90
34# if defined PIO_LIB && defined DISTRIBUTE
35 PRIVATE :: obs_read_pio
36# endif
37!
38 CONTAINS
39!
40!***********************************************************************
41 SUBROUTINE obs_read (ng, model, backward)
42!***********************************************************************
43!
44! Imported variable declarations.
45!
46 logical, intent(in) :: backward
47!
48 integer, intent(in) :: ng, model
49!
50! Local variable declarations.
51!
52 character (len=*), parameter :: myfile = &
53 & __FILE__
54!
55!-----------------------------------------------------------------------
56! Write out time-averaged fields according to IO type.
57!-----------------------------------------------------------------------
58!
59 SELECT CASE (obs(ng)%IOtype)
60 CASE (io_nf90)
61 CALL obs_read_nf90 (ng, model, backward)
62
63# if defined PIO_LIB && defined DISTRIBUTE
64 CASE (io_pio)
65 CALL obs_read_pio (ng, model, backward)
66# endif
67 CASE DEFAULT
68 IF (master) WRITE (stdout,10) obs(ng)%IOtype
69 exit_flag=2
70 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
71 END SELECT
72!
73 10 FORMAT (' OBS_READ - Illegal input file type, io_type = ',i0, &
74 & /,12x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
75!
76 RETURN
77 END SUBROUTINE obs_read
78!
79!***********************************************************************
80 SUBROUTINE obs_read_nf90 (ng, model, backward)
81!***********************************************************************
82!
83 USE mod_netcdf
84!
85! Imported variable declarations.
86!
87 logical, intent(in) :: backward
88!
89 integer, intent(in) :: ng, model
90!
91! Local variable declarations.
92!
93 logical :: readnlmod, readtlmod
94!
95 integer :: mstr, mend
96 integer :: i, iobs, itrc, status
97!
98 character (len=22) :: t_code
99
100 character (len=*), parameter :: myfile = &
101 & __FILE__//", obs_read_nf90"
102!
103 sourcefile=myfile
104!
105!---------------------------------------------------------------------
106! Read observation variables needed for interpolating the model
107! state at the observation locations.
108!---------------------------------------------------------------------
109!
110 IF (processobs(ng)) THEN
111# if defined TLM_OBS
112# if defined I4DVAR_ANA_SENSITIVITY
113 readnlmod=.false.
114# else
115 readnlmod=.true.
116# endif
117 readtlmod=.true.
118# else
119 readnlmod=.false.
120 readtlmod=.false.
121# endif
122!
123! Initialize observations processing counters.
124!
125 DO i=1,nobsvar(ng)
126 fourdvar(ng)%ObsCount(i)=0
127 fourdvar(ng)%ObsReject(i)=0
128 END DO
129 IF (backward) THEN
130 obssurvey(ng)=obssurvey(ng)-1
131 ELSE
132 obssurvey(ng)=obssurvey(ng)+1
133 END IF
134!
135! Set number of observations to process.
136!
137 nobs(ng)=fourdvar(ng)%NobsSurvey(obssurvey(ng))
138!
139! Set number of datum to process at current time-step.
140!
141 IF (backward) THEN
142 nendobs(ng)=nstrobs(ng)-1
143 nstrobs(ng)=nstrobs(ng)-nobs(ng)
144 ELSE
145 nstrobs(ng)=nendobs(ng)+1
146 nendobs(ng)=nstrobs(ng)+nobs(ng)-1
147 END IF
148!
149! Set starting index of obervation vectors for reading. In weak
150! constraint, the entire observation data is loaded. Otherwise,
151! only the observartion for the current time window are loaded
152! and started from vector index one.
153!
154# ifdef WEAK_CONSTRAINT
155 mstr=nstrobs(ng)
156 mend=nendobs(ng)
157# else
158 mstr=1
159 mend=nobs(ng)
160# endif
161!
162! Read in observation type identifier.
163!
164 CALL netcdf_get_ivar (ng, model, obs(ng)%name, vname(1,idotyp), &
165 & obstype(mstr:), &
166 & ncid = obs(ng)%ncid, &
167 & start = (/nstrobs(ng)/), &
168 & total = (/nobs(ng)/))
169 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
170!
171! Read in observation provenance identifier.
172!
173 CALL netcdf_get_ivar (ng, model, obs(ng)%name, vname(1,idopro), &
174 & obsprov(mstr:), &
175 & ncid = obs(ng)%ncid, &
176 & start = (/nstrobs(ng)/), &
177 & total = (/nobs(ng)/))
178 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
179!
180! Read in observation time (days).
181!
182 CALL netcdf_get_time (ng, model, obs(ng)%name, vname(1,idobst), &
183 & rclock%DateNumber, tobs(mstr:), &
184 & ncid = obs(ng)%ncid, &
185 & start = (/nstrobs(ng)/), &
186 & total = (/nobs(ng)/))
187 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
188!
189! Read in observation X-location (grid units).
190!
191 CALL netcdf_get_fvar (ng, model, obs(ng)%name, vname(1,idobsx), &
192 & xobs(mstr:), &
193 & ncid = obs(ng)%ncid, &
194 & start = (/nstrobs(ng)/), &
195 & total = (/nobs(ng)/))
196 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
197!
198! Read in observation Y-location (grid units).
199!
200 CALL netcdf_get_fvar (ng, model, obs(ng)%name, vname(1,idobsy), &
201 & yobs(mstr:), &
202 & ncid = obs(ng)%ncid, &
203 & start = (/nstrobs(ng)/), &
204 & total = (/nobs(ng)/))
205 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
206
207# ifdef SOLVE3D
208!
209! Read in observation depth, Zobs. If negative, depth is meter. If
210! greater than zero, Zobs is in model fractional vertical levels
211! (1 <= Zobs <= N). If Zobs < 0, its fractional level value is
212! computed in routine "extract_obs3d" and over-written so it can
213! be written into the observation NetCDF file for latter use.
214!
215 IF (wrote_zobs(ng)) THEN
216 CALL netcdf_get_fvar (ng, model, obs(ng)%name, &
217 & vname(1,idobsz), &
218 & zobs(mstr:), &
219 & ncid = obs(ng)%ncid, &
220 & start = (/nstrobs(ng)/), &
221 & total = (/nobs(ng)/))
222 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
223 ELSE
224 CALL netcdf_get_fvar (ng, model, obs(ng)%name, &
225 & vname(1,idobsd), &
226 & zobs(mstr:), &
227 & ncid = obs(ng)%ncid, &
228 & start = (/nstrobs(ng)/), &
229 & total = (/nobs(ng)/))
230 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
231 END IF
232 load_zobs(ng)=.false.
233 IF ((minval(zobs).lt.0.0_r8).or. &
234 & (maxval(zobs).lt.0.0_r8)) THEN
235 load_zobs(ng)=.true.
236 END IF
237
238# ifdef DISTRIBUTE
239!
240! If distributed-memory and SOME or ALL observations have Zobs in
241! meters (Zobs < 0), zero-out Zobs values in all nodes but itself
242! to facilitate exchages between tiles latter before writting into
243! observation NetCDF file.
244!
245 IF (.not.wrote_zobs(ng).and.load_zobs(ng)) THEN
246 CALL obs_depth (ng, myrank, model)
247 END IF
248# endif
249# endif
250!
251! Read in observation values.
252!
253 CALL netcdf_get_fvar (ng, model, obs(ng)%name, vname(1,idoval), &
254 & obsval(mstr:), &
255 & ncid = obs(ng)%ncid, &
256 & start = (/nstrobs(ng)/), &
257 & total = (/nobs(ng)/))
258 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
259!
260! If appropriate, read in observation meta values.
261!
262 IF (haveobsmeta(ng)) THEN
263 CALL netcdf_get_fvar (ng, model, obs(ng)%name, &
264 & vname(1,idomet), &
265 & obsmeta(mstr:), &
266 & ncid = obs(ng)%ncid, &
267 & start = (/nstrobs(ng)/), &
268 & total = (/nobs(ng)/))
269 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
270 END IF
271!
272# ifdef WEAK_CONSTRAINT
273! Read in observation error covariance.
274# else
275! Read in observation error covariance. To avoid successive divisions,
276! convert to inverse observation error covariance.
277# endif
278!
279 CALL netcdf_get_fvar (ng, model, obs(ng)%name, vname(1,idoerr), &
280 & obserr(mstr:), &
281 & ncid = obs(ng)%ncid, &
282 & start = (/nstrobs(ng)/), &
283 & total = (/nobs(ng)/))
284 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
285
286# ifndef WEAK_CONSTRAINT
287 DO iobs=1,nobs(ng)
288 obserr(iobs)=1.0_r8/obserr(iobs)
289 END DO
290# endif
291!
292! Read in nonlinear model values at observation locations.
293!
294 IF (readnlmod.and.havenlmod(ng)) THEN
295# if defined VERIFICATION || defined TLM_CHECK
296 CALL netcdf_get_fvar (ng, model, dav(ng)%name, &
297 & vname(1,idnlmo), &
298 & nlmodval(mstr:), &
299 & ncid = dav(ng)%ncid, &
300 & start = (/nstrobs(ng)/), &
301 & total = (/nobs(ng)/))
302 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
303# else
304 CALL netcdf_get_fvar (ng, model, dav(ng)%name, &
305 & vname(1,idnlmo), &
306 & nlmodval(mstr:), &
307 & ncid = dav(ng)%ncid, &
308 & start = (/nstrobs(ng),outer/), &
309 & total = (/nobs(ng),1/))
310 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
311# endif
312# ifdef BGQC_NOT_NEEDED
313 CALL netcdf_get_fvar (ng, model, dav(ng)%name, &
314 & vname(1,idbger), &
315 & bgerr(mstr:), &
316 & ncid = dav(ng)%ncid, &
317 & start = (/nstrobs(ng)/), &
318 & total = (/nobs(ng)/))
319 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
320 CALL netcdf_get_fvar (ng, model, dav(ng)%name, &
321 & vname(1,idobss), &
322 & obsscale(mstr:), &
323 & ncid = dav(ng)%ncid, &
324 & start = (/nstrobs(ng)/), &
325 & total = (/nobs(ng)/))
326 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
327# endif
328 END IF
329
330# if defined TLM_OBS && !defined WEAK_CONSTRAINT
331!
332! If adjoint pass and incremental 4DVar, read in tangent linear model
333! values at observation locations.
334!
335 IF (readtlmod.and.havetlmod(ng)) THEN
336 CALL netcdf_get_fvar (ng, model, dav(ng)%name, &
337 & vname(1,idtlmo), &
338 & tlmodval(mstr:), &
339 & ncid = dav(ng)%ncid, &
340 & start = (/nstrobs(ng)/), &
341 & total = (/nobs(ng)/))
342 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
343
344# if defined I4DVAR
345!
346! Reset TLM values to zero at the first pass of the inner loop.
347!
348 IF (inner.eq.0) THEN
349 DO iobs=1,mobs
350 tlmodval(iobs)=0.0_r8
351 END DO
352 END IF
353# endif
354 END IF
355# endif
356!
357!-----------------------------------------------------------------------
358! Set counters for number of observations to processed for each state
359! variable.
360!-----------------------------------------------------------------------
361!
362 DO iobs=mstr,mend
363 IF (obstype(iobs).eq.obsstate2type(isfsur)) THEN
364 fourdvar(ng)%ObsCount(isfsur)= &
365 & fourdvar(ng)%ObsCount(isfsur)+1
366 ELSE IF (obstype(iobs).eq.obsstate2type(isubar)) THEN
367 fourdvar(ng)%ObsCount(isubar)= &
368 & fourdvar(ng)%ObsCount(isubar)+1
369 ELSE IF (obstype(iobs).eq.obsstate2type(isvbar)) THEN
370 fourdvar(ng)%ObsCount(isvbar)= &
371 & fourdvar(ng)%ObsCount(isvbar)+1
372# ifdef SOLVE3D
373 ELSE IF (obstype(iobs).eq.obsstate2type(isuvel)) THEN
374 fourdvar(ng)%ObsCount(isuvel)= &
375 & fourdvar(ng)%ObsCount(isuvel)+1
376 ELSE IF (obstype(iobs).eq.obsstate2type(isvvel)) THEN
377 fourdvar(ng)%ObsCount(isvvel)= &
378 & fourdvar(ng)%ObsCount(isvvel)+1
379 ELSE IF (obstype(iobs).eq.obsstate2type(isradial)) THEN
380 haveobsmeta(ng)=.true.
381 fourdvar(ng)%ObsCount(isradial)= &
382 & fourdvar(ng)%ObsCount(isradial)+1
383 ELSE
384 DO itrc=1,nt(ng)
385 IF (obstype(iobs).eq.obsstate2type(istvar(itrc))) THEN
386 i=istvar(itrc)
387 fourdvar(ng)%ObsCount(i)=fourdvar(ng)%ObsCount(i)+1
388 END IF
389 END DO
390# endif
391 END IF
392 END DO
393!
394!-----------------------------------------------------------------------
395! If applicable, set next observation survey time to process.
396!-----------------------------------------------------------------------
397!
398 IF (master) THEN
399 CALL time_string (obstime(ng), t_code)
400 WRITE (stdout,10) obstime(ng)*sec2day, t_code
401 END IF
402 IF (backward) THEN
403 IF ((obssurvey(ng)-1).ge.1) THEN
404 obstime(ng)=fourdvar(ng)%SurveyTime(obssurvey(ng)-1)*day2sec
405 END IF
406 ELSE
407 IF ((obssurvey(ng)+1).le.nsurvey(ng)) THEN
408 obstime(ng)=fourdvar(ng)%SurveyTime(obssurvey(ng)+1)*day2sec
409 END IF
410 END IF
411 END IF
412!
413 10 FORMAT (/,' Number of State Observations Processed:',2x, &
414 & 'ObsTime = ',f12.4,',',t68,a,/,/, &
415 & 10x,'Variable',10x,'IstrObs',4x,'IendObs',6x,'Count', &
416 & 3x,'Rejected',/)
417!
418 RETURN
419 END SUBROUTINE obs_read_nf90
420
421# if defined PIO_LIB && defined DISTRIBUTE
422!
423!***********************************************************************
424 SUBROUTINE obs_read_pio (ng, model, backward)
425!***********************************************************************
426!
428!
429! Imported variable declarations.
430!
431 logical, intent(in) :: backward
432!
433 integer, intent(in) :: ng, model
434!
435! Local variable declarations.
436!
437 logical :: readnlmod, readtlmod
438!
439 integer :: mstr, mend
440 integer :: i, iobs, itrc, status
441!
442 character (len=22) :: t_code
443
444 character (len=*), parameter :: myfile = &
445 & __FILE__//", obs_read_pio"
446!
447 sourcefile=myfile
448!
449!---------------------------------------------------------------------
450! Read observation variables needed for interpolating the model
451! state at the observation locations.
452!---------------------------------------------------------------------
453!
454 IF (processobs(ng)) THEN
455# if defined TLM_OBS
456# if defined I4DVAR_ANA_SENSITIVITY
457 readnlmod=.false.
458# else
459 readnlmod=.true.
460# endif
461 readtlmod=.true.
462# else
463 readnlmod=.false.
464 readtlmod=.false.
465# endif
466!
467! Initialize observations processing counters.
468!
469 DO i=1,nobsvar(ng)
470 fourdvar(ng)%ObsCount(i)=0
471 fourdvar(ng)%ObsReject(i)=0
472 END DO
473 IF (backward) THEN
474 obssurvey(ng)=obssurvey(ng)-1
475 ELSE
476 obssurvey(ng)=obssurvey(ng)+1
477 END IF
478!
479! Set number of observations to process.
480!
481 nobs(ng)=fourdvar(ng)%NobsSurvey(obssurvey(ng))
482!
483! Set number of datum to process at current time-step.
484!
485 IF (backward) THEN
486 nendobs(ng)=nstrobs(ng)-1
487 nstrobs(ng)=nstrobs(ng)-nobs(ng)
488 ELSE
489 nstrobs(ng)=nendobs(ng)+1
490 nendobs(ng)=nstrobs(ng)+nobs(ng)-1
491 END IF
492!
493! Set starting index of obervation vectors for reading. In weak
494! constraint, the entire observation data is loaded. Otherwise,
495! only the observartion for the current time window are loaded
496! and started from vector index one.
497!
498# ifdef WEAK_CONSTRAINT
499 mstr=nstrobs(ng)
500 mend=nendobs(ng)
501# else
502 mstr=1
503 mend=nobs(ng)
504# endif
505!
506! Read in observation type identifier.
507!
508 CALL pio_netcdf_get_ivar (ng, model, obs(ng)%name, &
509 & vname(1,idotyp), &
510 & obstype(mstr:), &
511 & piofile = obs(ng)%pioFile, &
512 & start = (/nstrobs(ng)/), &
513 & total = (/nobs(ng)/))
514 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
515!
516! Read in observation provenance identifier.
517!
518 CALL pio_netcdf_get_ivar (ng, model, obs(ng)%name, &
519 & vname(1,idopro), &
520 & obsprov(mstr:), &
521 & piofile = obs(ng)%pioFile, &
522 & start = (/nstrobs(ng)/), &
523 & total = (/nobs(ng)/))
524 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
525!
526! Read in observation time (days).
527!
528 CALL pio_netcdf_get_time (ng, model, obs(ng)%name, &
529 & vname(1,idobst), &
530 & rclock%DateNumber, tobs(mstr:), &
531 & piofile = obs(ng)%pioFile, &
532 & start = (/nstrobs(ng)/), &
533 & total = (/nobs(ng)/))
534 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
535!
536! Read in observation X-location (grid units).
537!
538 CALL pio_netcdf_get_fvar (ng, model, obs(ng)%name, &
539 & vname(1,idobsx), &
540 & xobs(mstr:), &
541 & piofile = obs(ng)%pioFile, &
542 & start = (/nstrobs(ng)/), &
543 & total = (/nobs(ng)/))
544 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
545!
546! Read in observation Y-location (grid units).
547!
548 CALL pio_netcdf_get_fvar (ng, model, obs(ng)%name, &
549 & vname(1,idobsy), &
550 & yobs(mstr:), &
551 & piofile = obs(ng)%pioFile, &
552 & start = (/nstrobs(ng)/), &
553 & total = (/nobs(ng)/))
554 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
555
556# ifdef SOLVE3D
557!
558! Read in observation depth, Zobs. If negative, depth is meter. If
559! greater than zero, Zobs is in model fractional vertical levels
560! (1 <= Zobs <= N). If Zobs < 0, its fractional level value is
561! computed in routine "extract_obs3d" and over-written so it can
562! be written into the observation NetCDF file for latter use.
563!
564 IF (wrote_zobs(ng)) THEN
565 CALL pio_netcdf_get_fvar (ng, model, obs(ng)%name, &
566 & vname(1,idobsz), &
567 & zobs(mstr:), &
568 & piofile = obs(ng)%pioFile, &
569 & start = (/nstrobs(ng)/), &
570 & total = (/nobs(ng)/))
571 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
572 ELSE
573 CALL pio_netcdf_get_fvar (ng, model, obs(ng)%name, &
574 & vname(1,idobsd), &
575 & zobs(mstr:), &
576 & piofile = obs(ng)%pioFile, &
577 & start = (/nstrobs(ng)/), &
578 & total = (/nobs(ng)/))
579 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
580 END IF
581 load_zobs(ng)=.false.
582 IF ((minval(zobs).lt.0.0_r8).or. &
583 & (maxval(zobs).lt.0.0_r8)) THEN
584 load_zobs(ng)=.true.
585 END IF
586!
587! If distributed-memory and SOME or ALL observations have Zobs in
588! meters (Zobs < 0), zero-out Zobs values in all nodes but itself
589! to facilitate exchages between tiles latter before writting into
590! observation NetCDF file.
591!
592 IF (.not.wrote_zobs(ng).and.load_zobs(ng)) THEN
593 CALL obs_depth (ng, myrank, model)
594 END IF
595# endif
596!
597! Read in observation values.
598!
599 CALL pio_netcdf_get_fvar (ng, model, obs(ng)%name, &
600 & vname(1,idoval), &
601 & obsval(mstr:), &
602 & piofile = obs(ng)%pioFile, &
603 & start = (/nstrobs(ng)/), &
604 & total = (/nobs(ng)/))
605 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
606!
607! If appropriate, read in observation meta values.
608!
609 IF (haveobsmeta(ng)) THEN
610 CALL pio_netcdf_get_fvar (ng, model, obs(ng)%name, &
611 & vname(1,idomet), &
612 & obsmeta(mstr:), &
613 & piofile = obs(ng)%pioFile, &
614 & start = (/nstrobs(ng)/), &
615 & total = (/nobs(ng)/))
616 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
617 END IF
618!
619# ifdef WEAK_CONSTRAINT
620! Read in observation error covariance.
621# else
622! Read in observation error covariance. To avoid successive divisions,
623! convert to inverse observation error covariance.
624# endif
625!
626 CALL pio_netcdf_get_fvar (ng, model, obs(ng)%name, &
627 & vname(1,idoerr), &
628 & obserr(mstr:), &
629 & piofile = obs(ng)%pioFile, &
630 & start = (/nstrobs(ng)/), &
631 & total = (/nobs(ng)/))
632 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
633
634# ifndef WEAK_CONSTRAINT
635 DO iobs=1,nobs(ng)
636 obserr(iobs)=1.0_r8/obserr(iobs)
637 END DO
638# endif
639!
640! Read in nonlinear model values at observation locations.
641!
642 IF (readnlmod.and.havenlmod(ng)) THEN
643# ifdef VERIFICATION
644 CALL pio_netcdf_get_fvar (ng, model, dav(ng)%name, &
645 & vname(1,idnlmo), &
646 & nlmodval(mstr:), &
647 & piofile = dav(ng)%pioFile, &
648 & start = (/nstrobs(ng)/), &
649 & total = (/nobs(ng)/))
650 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
651# else
652 CALL pio_netcdf_get_fvar (ng, model, dav(ng)%name, &
653 & vname(1,idnlmo), &
654 & nlmodval(mstr:), &
655 & piofile = dav(ng)%pioFile, &
656 & start = (/nstrobs(ng),outer/), &
657 & total = (/nobs(ng),1/))
658 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
659# endif
660# ifdef BGQC_NOT_NEEDED
661 CALL pio_netcdf_get_fvar (ng, model, dav(ng)%name, &
662 & vname(1,idbger), &
663 & bgerr(mstr:), &
664 & piofile = dav(ng)%pioFile, &
665 & start = (/nstrobs(ng)/), &
666 & total = (/nobs(ng)/))
667 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
668
669 CALL pio_netcdf_get_fvar (ng, model, dav(ng)%name, &
670 & vname(1,idobss), &
671 & obsscale(mstr:), &
672 & piofile = dav(ng)%pioFile, &
673 & start = (/nstrobs(ng)/), &
674 & total = (/nobs(ng)/))
675 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
676# endif
677 END IF
678
679# if defined TLM_OBS && !defined WEAK_CONSTRAINT
680!
681! If adjoint pass and incremental 4DVar, read in tangent linear model
682! values at observation locations.
683!
684 IF (readtlmod.and.havetlmod(ng)) THEN
685 CALL pio_netcdf_get_fvar (ng, model, dav(ng)%name, &
686 & vname(1,idtlmo), &
687 & tlmodval(mstr:), &
688 & piofile = dav(ng)%pioFile, &
689 & start = (/nstrobs(ng)/), &
690 & total = (/nobs(ng)/))
691 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
692
693# if defined I4DVAR
694!
695! Reset TLM values to zero at the first pass of the inner loop.
696!
697 IF (inner.eq.0) THEN
698 DO iobs=1,mobs
699 tlmodval(iobs)=0.0_r8
700 END DO
701 END IF
702# endif
703 END IF
704# endif
705!
706!-----------------------------------------------------------------------
707! Set counters for number of observations to processed for each state
708! variable.
709!-----------------------------------------------------------------------
710!
711 DO iobs=mstr,mend
712 IF (obstype(iobs).eq.obsstate2type(isfsur)) THEN
713 fourdvar(ng)%ObsCount(isfsur)= &
714 & fourdvar(ng)%ObsCount(isfsur)+1
715 ELSE IF (obstype(iobs).eq.obsstate2type(isubar)) THEN
716 fourdvar(ng)%ObsCount(isubar)= &
717 & fourdvar(ng)%ObsCount(isubar)+1
718 ELSE IF (obstype(iobs).eq.obsstate2type(isvbar)) THEN
719 fourdvar(ng)%ObsCount(isvbar)= &
720 & fourdvar(ng)%ObsCount(isvbar)+1
721# ifdef SOLVE3D
722 ELSE IF (obstype(iobs).eq.obsstate2type(isuvel)) THEN
723 fourdvar(ng)%ObsCount(isuvel)= &
724 & fourdvar(ng)%ObsCount(isuvel)+1
725 ELSE IF (obstype(iobs).eq.obsstate2type(isvvel)) THEN
726 fourdvar(ng)%ObsCount(isvvel)= &
727 & fourdvar(ng)%ObsCount(isvvel)+1
728 ELSE IF (obstype(iobs).eq.obsstate2type(isradial)) THEN
729 haveobsmeta(ng)=.true.
730 fourdvar(ng)%ObsCount(isradial)= &
731 & fourdvar(ng)%ObsCount(isradial)+1
732 ELSE
733 DO itrc=1,nt(ng)
734 IF (obstype(iobs).eq.obsstate2type(istvar(itrc))) THEN
735 i=istvar(itrc)
736 fourdvar(ng)%ObsCount(i)=fourdvar(ng)%ObsCount(i)+1
737 END IF
738 END DO
739# endif
740 END IF
741 END DO
742!
743!-----------------------------------------------------------------------
744! If applicable, set next observation survey time to process.
745!-----------------------------------------------------------------------
746!
747 IF (master) THEN
748 CALL time_string (obstime(ng), t_code)
749 WRITE (stdout,10) obstime(ng)*sec2day, t_code
750 END IF
751 IF (backward) THEN
752 IF ((obssurvey(ng)-1).ge.1) THEN
753 obstime(ng)=fourdvar(ng)%SurveyTime(obssurvey(ng)-1)*day2sec
754 END IF
755 ELSE
756 IF ((obssurvey(ng)+1).le.nsurvey(ng)) THEN
757 obstime(ng)=fourdvar(ng)%SurveyTime(obssurvey(ng)+1)*day2sec
758 END IF
759 END IF
760 END IF
761!
762 10 FORMAT (/,' Number of State Observations Processed:',2x, &
763 & 'ObsTime = ',f12.4,',',t68,a,/,/, &
764 & 10x,'Variable',10x,'IstrObs',4x,'IendObs',6x,'Count', &
765 & 3x,'Rejected',/)
766!
767 RETURN
768 END SUBROUTINE obs_read_pio
769# endif
770#endif
771 END MODULE obs_read_mod
subroutine, public time_string(mytime, date_string)
Definition dateclock.F:1272
type(t_fourdvar), dimension(:), allocatable fourdvar
integer, dimension(:), allocatable obsprov
integer, dimension(:), allocatable nobs
integer, dimension(:), allocatable nobsvar
real(r8), dimension(:), allocatable obsval
logical, dimension(:), allocatable havenlmod
real(r8), dimension(:), allocatable bgerr
logical, dimension(:), allocatable wrote_zobs
real(r8), dimension(:), allocatable obsscale
real(r8), dimension(:), allocatable obserr
real(r8), dimension(:), allocatable obsmeta
logical, dimension(:), allocatable load_zobs
integer, dimension(:), allocatable obstype
logical, dimension(:), allocatable haveobsmeta
logical, dimension(:), allocatable havetlmod
integer, dimension(:), allocatable nsurvey
real(r8), dimension(:), allocatable zobs
real(r8), dimension(:), allocatable nlmodval
real(dp), dimension(:), allocatable tobs
integer, dimension(:), allocatable obssurvey
real(r8), dimension(:), allocatable xobs
integer, dimension(:), allocatable obsstate2type
real(r8), dimension(:), allocatable yobs
integer, dimension(:), allocatable nstrobs
logical, dimension(:), allocatable processobs
integer, dimension(:), allocatable nendobs
type(t_io), dimension(:), allocatable obs
type(t_io), dimension(:), allocatable dav
integer stdout
character(len=256) sourcefile
integer idoerr
integer idnlmo
integer idobss
integer idtlmo
integer idobst
integer, parameter io_nf90
Definition mod_ncparam.F:95
integer idopro
integer idobsy
integer idobsd
integer isvvel
integer, parameter io_pio
Definition mod_ncparam.F:96
integer isvbar
integer idomet
integer idoval
integer, dimension(:), allocatable istvar
integer idobsz
integer isuvel
integer isfsur
character(len=maxlen), dimension(6, 0:nv) vname
integer isubar
integer idbger
integer idobsx
integer isradial
integer idotyp
logical master
integer, dimension(:), allocatable nt
Definition mod_param.F:489
real(dp), parameter day2sec
real(dp), dimension(:), allocatable obstime
type(t_clock) rclock
real(dp), parameter sec2day
integer exit_flag
integer inner
integer noerror
integer outer
subroutine, private obs_read_pio(ng, model, backward)
Definition obs_read.F:425
subroutine, private obs_read_nf90(ng, model, backward)
Definition obs_read.F:81
subroutine, public obs_read(ng, model, backward)
Definition obs_read.F:42
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine obs_depth(ng, tile, model)
Definition obs_depth.F:6