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

Functions/Subroutines

subroutine, public obs_initial (ng, model, backward)
 
subroutine, private obs_initial_nf90 (ng, model, backward)
 
subroutine, private obs_initial_pio (ng, model, backward)
 

Function/Subroutine Documentation

◆ obs_initial()

subroutine, public obs_initial_mod::obs_initial ( integer, intent(in) ng,
integer, intent(in) model,
logical, intent(in) backward )

Definition at line 40 of file obs_initial.F.

41!***********************************************************************
42!
43! Imported variable declarations.
44!
45 logical, intent(in) :: backward
46!
47 integer, intent(in) :: ng, model
48!
49! Local variable declarations.
50!
51 character (len=*), parameter :: MyFile = &
52 & __FILE__
53!
54!-----------------------------------------------------------------------
55! Write out time-averaged fields according to IO type.
56!-----------------------------------------------------------------------
57!
58 SELECT CASE (obs(ng)%IOtype)
59 CASE (io_nf90)
60 CALL obs_initial_nf90 (ng, model, backward)
61
62# if defined PIO_LIB && defined DISTRIBUTE
63 CASE (io_pio)
64 CALL obs_initial_pio (ng, model, backward)
65# endif
66 CASE DEFAULT
67 IF (master) WRITE (stdout,10) obs(ng)%IOtype
68 exit_flag=2
69 END SELECT
70 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
71!
72 10 FORMAT (' OBS_INITIAL - Illegal input file type, io_type = ',i0, &
73 & /,15x,'Check KeyWord ''INP_LIB'' in ''roms.in''.')
74!
75 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::obs, obs_initial_nf90(), obs_initial_pio(), and mod_iounits::stdout.

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

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

◆ obs_initial_nf90()

subroutine, private obs_initial_mod::obs_initial_nf90 ( integer, intent(in) ng,
integer, intent(in) model,
logical, intent(in) backward )
private

Definition at line 79 of file obs_initial.F.

80!***********************************************************************
81!
82 USE mod_netcdf
83!
84! Imported variable declarations.
85!
86 logical, intent(in) :: backward
87!
88 integer, intent(in) :: ng, model
89!
90! Local variable declarations.
91!
92 logical, dimension(NV) :: got_var(NV)
93!
94 integer :: Ifirst, i, nvd, recdim, status
95 integer :: Vsize(4)
96!
97# ifdef WEAK_CONSTRAINT
98 real(r8), parameter :: IniVal = 0.0_r8
99# endif
100 real(r8) :: tend
101!
102 character (len=*), parameter :: MyFile = &
103 & __FILE__//"obs_initial_nf90"
104!
105 sourcefile=myfile
106!
107!-----------------------------------------------------------------------
108! Inquire about the contents of observation NetCDF file: Inquire about
109! the dimensions and variables.
110!-----------------------------------------------------------------------
111!
112 query : IF (obs(ng)%ncid.eq.-1) THEN
113!
114! Open observations NetCDF file.
115!
116 CALL netcdf_open (ng, model, obs(ng)%name, 1, obs(ng)%ncid)
117 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
118 WRITE (stdout,10) trim(obs(ng)%name)
119 RETURN
120 END IF
121!
122! Inquire about the variables.
123!
124 CALL netcdf_inq_var (ng, model, obs(ng)%name, &
125 & obs(ng)%ncid)
126 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
127!
128! Initialize logical switches.
129!
130 DO i=1,nv
131 got_var(i)=.false.
132 END DO
133!
134! Scan variable list from observation NetCDF and activate switches for
135! required variables.
136!
137 DO i=1,n_var
138 IF (trim(var_name(i)).eq.trim(vname(1,idoday))) THEN
139 got_var(idoday)=.true.
140 obs(ng)%Vid(idoday)=var_id(i)
141 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idnobs))) THEN
142 got_var(idnobs)=.true.
143 obs(ng)%Vid(idnobs)=var_id(i)
144 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idotyp))) THEN
145 got_var(idotyp)=.true.
146 obs(ng)%Vid(idotyp)=var_id(i)
147 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobst))) THEN
148 got_var(idobst)=.true.
149 obs(ng)%Vid(idobst)=var_id(i)
150# ifdef SOLVE3D
151 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsd))) THEN
152 got_var(idobsd)=.true.
153 obs(ng)%Vid(idobsd)=var_id(i)
154# endif
155 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsx))) THEN
156 got_var(idobsx)=.true.
157 obs(ng)%Vid(idobsx)=var_id(i)
158 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsy))) THEN
159 got_var(idobsy)=.true.
160 obs(ng)%Vid(idobsy)=var_id(i)
161# ifdef SOLVE3D
162 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsz))) THEN
163 got_var(idobsz)=.true.
164 obs(ng)%Vid(idobsz)=var_id(i)
165# endif
166 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idoerr))) THEN
167 got_var(idoerr)=.true.
168 obs(ng)%Vid(idoerr)=var_id(i)
169 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idoval))) THEN
170 got_var(idoval)=.true.
171 obs(ng)%Vid(idoval)=var_id(i)
172 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idomet))) THEN
173 got_var(idomet)=.true.
174 haveobsmeta(ng)=.true.
175 obs(ng)%Vid(idomet)=var_id(i)
176 END IF
177 END DO
178!
179! Check if needed obsrvation variables are available.
180!
181 IF (.not.got_var(idoday)) THEN
182 IF (master) WRITE (stdout,20) trim(vname(1,idoday)), &
183 & trim(obs(ng)%name)
184 exit_flag=2
185 RETURN
186 END IF
187 IF (.not.got_var(idnobs)) THEN
188 IF (master) WRITE (stdout,20) trim(vname(1,idnobs)), &
189 & trim(obs(ng)%name)
190 exit_flag=2
191 RETURN
192 END IF
193 IF (.not.got_var(idotyp)) THEN
194 IF (master) WRITE (stdout,20) trim(vname(1,idotyp)), &
195 & trim(obs(ng)%name)
196 exit_flag=2
197 RETURN
198 END IF
199 IF (.not.got_var(idobst)) THEN
200 IF (master) WRITE (stdout,20) trim(vname(1,idobst)), &
201 & trim(obs(ng)%name)
202 exit_flag=2
203 RETURN
204 END IF
205# ifdef SOLVE3D
206 IF (.not.got_var(idobsd)) THEN
207 IF (master) WRITE (stdout,20) trim(vname(1,idobsd)), &
208 & trim(obs(ng)%name)
209 exit_flag=2
210 RETURN
211 END IF
212# endif
213 IF (.not.got_var(idobsx)) THEN
214 IF (master) WRITE (stdout,20) trim(vname(1,idobsx)), &
215 & trim(obs(ng)%name)
216 exit_flag=2
217 RETURN
218 END IF
219 IF (.not.got_var(idobsy)) THEN
220 IF (master) WRITE (stdout,20) trim(vname(1,idobsy)), &
221 & trim(obs(ng)%name)
222 exit_flag=2
223 RETURN
224 END IF
225# ifdef SOLVE3D
226 IF (.not.got_var(idobsz)) THEN
227 IF (master) WRITE (stdout,20) trim(vname(1,idobsz)), &
228 & trim(obs(ng)%name)
229 exit_flag=2
230 RETURN
231 END IF
232# endif
233 IF (.not.got_var(idoerr)) THEN
234 IF (master) WRITE (stdout,20) trim(vname(1,idoerr)), &
235 & trim(obs(ng)%name)
236 exit_flag=2
237 RETURN
238 END IF
239 IF (.not.got_var(idoval)) THEN
240 IF (master) WRITE (stdout,20) trim(vname(1,idoval)), &
241 & trim(obs(ng)%name)
242 exit_flag=2
243 RETURN
244 END IF
245 END IF query
246!
247!-----------------------------------------------------------------------
248! Set observation processing variables.
249!-----------------------------------------------------------------------
250!
251! Determine if there is any data available in the model time
252! window. Set first survey record to process.
253!
254 IF (backward) THEN
255 ifirst=0
256# ifdef GENERIC_DSTART
257 tend=(time(ng)-ntimes(ng)*dt(ng))*sec2day
258# else
259 tend=(time(ng)-(ntstart(ng)-1)*dt(ng))*sec2day
260# endif
261 DO i=1,nsurvey(ng)
262 IF ((tend.le.fourdvar(ng)%SurveyTime(i)).and. &
263 & (fourdvar(ng)%SurveyTime(i).le.tdays(ng))) THEN
264 ifirst=max(ifirst,i)
265 END IF
266 END DO
267# ifndef SP4DVAR
268 IF (ifirst.eq.0) THEN
269 IF (master) WRITE (stdout,30) tend, tdays(ng)
270 exit_flag=2
271 RETURN
272 END IF
273# endif
274 ELSE
275 ifirst=nsurvey(ng)
276 tend=(time(ng)+ntimes(ng)*dt(ng))*sec2day
277 DO i=1,nsurvey(ng)
278 IF ((tdays(ng).le.fourdvar(ng)%SurveyTime(i)).and. &
279 & (fourdvar(ng)%SurveyTime(i).le.tend)) THEN
280 ifirst=min(ifirst,i)
281 END IF
282 END DO
283# ifndef SP4DVAR
284 IF (ifirst.eq.0) THEN
285 IF (master) WRITE (stdout,30) tdays(ng), tend
286 exit_flag=2
287 RETURN
288 END IF
289# endif
290 END IF
291 obstime(ng)=fourdvar(ng)%SurveyTime(ifirst)*day2sec
292!
293! Initialize observation survey counter. This is the counter of data
294! assimilation cycles or observations survey times on which the model
295! state is extracted (interpolated) at the observation locations.
296!
297 IF (backward) THEN
298 obssurvey(ng)=ifirst+1
299 ELSE
300 obssurvey(ng)=ifirst-1
301 END IF
302!
303! Initialize time switch to process model state at observation
304! locations.
305!
306 processobs(ng)=.false.
307
308# ifdef I4DVAR
309!
310! Initialize cost function misfit between model and observations.
311! The IF statement is to avoid rewritting its value before it is
312! written into the initial NetCDF file.
313!
314 IF (.not.backward) THEN
315 DO i=0,nobsvar(ng)
316 fourdvar(ng)%ObsCost(i)=0.0_r8
317 END DO
318 END IF
319# endif
320!
321! Set staring and ending observation indices.
322!
323 IF (backward) THEN
324 nstrobs(ng)=0
325 nendobs(ng)=0
326 DO i=1,ifirst
327 nstrobs(ng)=nstrobs(ng)+fourdvar(ng)%NobsSurvey(i)
328 END DO
329 nstrobs(ng)=nstrobs(ng)+1
330 ELSE
331 IF (ifirst.eq.1) THEN
332 nstrobs(ng)=0
333 nendobs(ng)=0
334 ELSE
335 nstrobs(ng)=0
336 nendobs(ng)=0
337 DO i=1,ifirst-1
338 nendobs(ng)=nendobs(ng)+fourdvar(ng)%NobsSurvey(i)
339 END DO
340 END IF
341 END IF
342!
343! Initialize total observation counter in structure array. Notice
344! that the zero index carries the total summation.
345!
346 fourdvar(ng)%ObsCount(0)=0
347 fourdvar(ng)%ObsReject(0)=0
348
349# ifdef WEAK_CONSTRAINT
350!
351! Initialize model values at observation locations.
352!
353 DO i=1,mobs
354 nlmodval(i)=inival
355 obsvetting(i)=inival
356# ifndef SP4DVAR
357 tlmodval(i)=inival
358# endif
359# if defined SP4DVAR && defined DISJOINTED && \
360 defined concurrent_kernel
361 tlmodval(i)=inival
362# endif
363 END DO
364# if defined SP4DVAR && defined DISJOINTED && \
365 !defined CONCURRENT_KERNEL
366 IF (model.ne.iadm) THEN
367 DO i=1,mobs
368 tlmodval(i)=inival
369 END DO
370 END IF
371# endif
372# endif
373!
374 10 FORMAT (/,' OBS_INITIAL_NF90 - unable to open input NetCDF', &
375 & ' file: ',a)
376 20 FORMAT (/,' OBS_INITIAL_NF90 - unable to find model variable:', &
377 & 1x,a,/,20x,'in input NetCDF file: ',a)
378 30 FORMAT (/,' OBS_INITIAL_NF90 - No are observations available', &
379 & ' for time window (days): ',/,12x,f12.4,' - ',f12.4,/)
380!
381 RETURN
subroutine, public netcdf_open(ng, model, ncname, omode, ncid)
character(len=100), dimension(mvars) var_name
Definition mod_netcdf.F:169
integer, dimension(mvars) var_id
Definition mod_netcdf.F:160
integer n_var
Definition mod_netcdf.F:152
subroutine, public netcdf_inq_var(ng, model, ncname, ncid, myvarname, searchvar, varid, nvardim, nvaratt)

References mod_scalars::day2sec, mod_scalars::dt, mod_scalars::exit_flag, strings_mod::founderror(), mod_fourdvar::fourdvar, mod_fourdvar::haveobsmeta, mod_param::iadm, mod_ncparam::idnobs, mod_ncparam::idobsd, mod_ncparam::idobst, mod_ncparam::idobsx, mod_ncparam::idobsy, mod_ncparam::idobsz, mod_ncparam::idoday, mod_ncparam::idoerr, mod_ncparam::idomet, mod_ncparam::idotyp, mod_ncparam::idoval, mod_parallel::master, mod_fourdvar::mobs, mod_netcdf::n_var, mod_fourdvar::nendobs, mod_netcdf::netcdf_inq_var(), mod_netcdf::netcdf_open(), mod_fourdvar::nlmodval, mod_fourdvar::nobsvar, mod_scalars::noerror, mod_fourdvar::nstrobs, mod_fourdvar::nsurvey, mod_scalars::ntimes, mod_scalars::ntstart, mod_ncparam::nv, mod_iounits::obs, mod_fourdvar::obssurvey, mod_scalars::obstime, mod_fourdvar::obsvetting, mod_fourdvar::processobs, mod_scalars::sec2day, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::tdays, mod_scalars::time, mod_netcdf::var_id, mod_netcdf::var_name, and mod_ncparam::vname.

Referenced by obs_initial().

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

◆ obs_initial_pio()

subroutine, private obs_initial_mod::obs_initial_pio ( integer, intent(in) ng,
integer, intent(in) model,
logical, intent(in) backward )
private

Definition at line 387 of file obs_initial.F.

388!***********************************************************************
389!
391!
392! Imported variable declarations.
393!
394 logical, intent(in) :: backward
395!
396 integer, intent(in) :: ng, model
397!
398! Local variable declarations.
399!
400 logical, dimension(NV) :: got_var(NV)
401!
402 integer :: Ifirst, i, nvd, recdim, status
403 integer :: Vsize(4)
404!
405# ifdef WEAK_CONSTRAINT
406 real(r8), parameter :: IniVal = 0.0_r8
407# endif
408 real(r8) :: tend
409!
410 character (len=*), parameter :: MyFile = &
411 & __FILE__//"obs_initial_nf90"
412!
413 sourcefile=myfile
414!
415!-----------------------------------------------------------------------
416! Inquire about the contents of observation NetCDF file: Inquire about
417! the dimensions and variables.
418!-----------------------------------------------------------------------
419!
420 query : IF (obs(ng)%ncid.eq.-1) THEN
421!
422! Open observations NetCDF file.
423!
424 CALL pio_netcdf_open (ng, model, obs(ng)%name, 1, &
425 & obs(ng)%pioFile)
426 IF (founderror(exit_flag, noerror, __line__, myfile)) THEN
427 WRITE (stdout,10) trim(obs(ng)%name)
428 RETURN
429 END IF
430!
431! Inquire about the variables.
432!
433 CALL pio_netcdf_inq_var (ng, model, obs(ng)%name, &
434 & obs(ng)%pioFile)
435 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
436!
437! Initialize logical switches.
438!
439 DO i=1,nv
440 got_var(i)=.false.
441 END DO
442!
443! Scan variable list from observation NetCDF and activate switches for
444! required variables.
445!
446 DO i=1,n_var
447 IF (trim(var_name(i)).eq.trim(vname(1,idoday))) THEN
448 got_var(idoday)=.true.
449 obs(ng)%pioVar(idoday)%vd=var_desc(i)
450 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idnobs))) THEN
451 got_var(idnobs)=.true.
452 obs(ng)%pioVar(idnobs)%vd=var_desc(i)
453 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idotyp))) THEN
454 got_var(idotyp)=.true.
455 obs(ng)%pioVar(idotyp)%vd=var_desc(i)
456 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobst))) THEN
457 got_var(idobst)=.true.
458 obs(ng)%pioVar(idobst)%vd=var_desc(i)
459# ifdef SOLVE3D
460 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsd))) THEN
461 got_var(idobsd)=.true.
462 obs(ng)%pioVar(idobsd)%vd=var_desc(i)
463# endif
464 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsx))) THEN
465 got_var(idobsx)=.true.
466 obs(ng)%pioVar(idobsx)%vd=var_desc(i)
467 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsy))) THEN
468 got_var(idobsy)=.true.
469 obs(ng)%pioVar(idobsy)%vd=var_desc(i)
470# ifdef SOLVE3D
471 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idobsz))) THEN
472 got_var(idobsz)=.true.
473 obs(ng)%pioVar(idobsz)%vd=var_desc(i)
474# endif
475 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idoerr))) THEN
476 got_var(idoerr)=.true.
477 obs(ng)%pioVar(idoerr)%vd=var_desc(i)
478 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idoval))) THEN
479 got_var(idoval)=.true.
480 obs(ng)%pioVar(idoval)%vd=var_desc(i)
481 ELSE IF (trim(var_name(i)).eq.trim(vname(1,idomet))) THEN
482 got_var(idomet)=.true.
483 haveobsmeta(ng)=.true.
484 obs(ng)%pioVar(idomet)%vd=var_desc(i)
485 END IF
486 END DO
487!
488! Check if needed obsrvation variables are available.
489!
490 IF (.not.got_var(idoday)) THEN
491 IF (master) WRITE (stdout,20) trim(vname(1,idoday)), &
492 & trim(obs(ng)%name)
493 exit_flag=2
494 RETURN
495 END IF
496 IF (.not.got_var(idnobs)) THEN
497 IF (master) WRITE (stdout,20) trim(vname(1,idnobs)), &
498 & trim(obs(ng)%name)
499 exit_flag=2
500 RETURN
501 END IF
502 IF (.not.got_var(idotyp)) THEN
503 IF (master) WRITE (stdout,20) trim(vname(1,idotyp)), &
504 & trim(obs(ng)%name)
505 exit_flag=2
506 RETURN
507 END IF
508 IF (.not.got_var(idobst)) THEN
509 IF (master) WRITE (stdout,20) trim(vname(1,idobst)), &
510 & trim(obs(ng)%name)
511 exit_flag=2
512 RETURN
513 END IF
514# ifdef SOLVE3D
515 IF (.not.got_var(idobsd)) THEN
516 IF (master) WRITE (stdout,20) trim(vname(1,idobsd)), &
517 & trim(obs(ng)%name)
518 exit_flag=2
519 RETURN
520 END IF
521# endif
522 IF (.not.got_var(idobsx)) THEN
523 IF (master) WRITE (stdout,20) trim(vname(1,idobsx)), &
524 & trim(obs(ng)%name)
525 exit_flag=2
526 RETURN
527 END IF
528 IF (.not.got_var(idobsy)) THEN
529 IF (master) WRITE (stdout,20) trim(vname(1,idobsy)), &
530 & trim(obs(ng)%name)
531 exit_flag=2
532 RETURN
533 END IF
534# ifdef SOLVE3D
535 IF (.not.got_var(idobsz)) THEN
536 IF (master) WRITE (stdout,20) trim(vname(1,idobsz)), &
537 & trim(obs(ng)%name)
538 exit_flag=2
539 RETURN
540 END IF
541# endif
542 IF (.not.got_var(idoerr)) THEN
543 IF (master) WRITE (stdout,20) trim(vname(1,idoerr)), &
544 & trim(obs(ng)%name)
545 exit_flag=2
546 RETURN
547 END IF
548 IF (.not.got_var(idoval)) THEN
549 IF (master) WRITE (stdout,20) trim(vname(1,idoval)), &
550 & trim(obs(ng)%name)
551 exit_flag=2
552 RETURN
553 END IF
554 END IF query
555!
556!-----------------------------------------------------------------------
557! Set observation processing variables.
558!-----------------------------------------------------------------------
559!
560! Determine if there is any data available in the model time
561! window. Set first survey record to process.
562!
563 IF (backward) THEN
564 ifirst=0
565# ifdef GENERIC_DSTART
566 tend=(time(ng)-ntimes(ng)*dt(ng))*sec2day
567# else
568 tend=(time(ng)-(ntstart(ng)-1)*dt(ng))*sec2day
569# endif
570 DO i=1,nsurvey(ng)
571 IF ((tend.le.fourdvar(ng)%SurveyTime(i)).and. &
572 & (fourdvar(ng)%SurveyTime(i).le.tdays(ng))) THEN
573 ifirst=max(ifirst,i)
574 END IF
575 END DO
576# ifndef SP4DVAR
577 IF (ifirst.eq.0) THEN
578 IF (master) WRITE (stdout,30) tend, tdays(ng)
579 exit_flag=2
580 RETURN
581 END IF
582# endif
583 ELSE
584 ifirst=nsurvey(ng)
585 tend=(time(ng)+ntimes(ng)*dt(ng))*sec2day
586 DO i=1,nsurvey(ng)
587 IF ((tdays(ng).le.fourdvar(ng)%SurveyTime(i)).and. &
588 & (fourdvar(ng)%SurveyTime(i).le.tend)) THEN
589 ifirst=min(ifirst,i)
590 END IF
591 END DO
592# ifndef SP4DVAR
593 IF (ifirst.eq.0) THEN
594 IF (master) WRITE (stdout,30) tdays(ng), tend
595 exit_flag=2
596 RETURN
597 END IF
598# endif
599 END IF
600 obstime(ng)=fourdvar(ng)%SurveyTime(ifirst)*day2sec
601!
602! Initialize observation survey counter. This is the counter of data
603! assimilation cycles or observations survey times on which the model
604! state is extracted (interpolated) at the observation locations.
605!
606 IF (backward) THEN
607 obssurvey(ng)=ifirst+1
608 ELSE
609 obssurvey(ng)=ifirst-1
610 END IF
611!
612! Initialize time switch to process model state at observation
613! locations.
614!
615 processobs(ng)=.false.
616
617# ifdef I4DVAR
618!
619! Initialize cost function misfit between model and observations.
620! The IF statement is to avoid rewritting its value before it is
621! written into the initial NetCDF file.
622!
623 IF (.not.backward) THEN
624 DO i=0,nobsvar(ng)
625 fourdvar(ng)%ObsCost(i)=0.0_r8
626 END DO
627 END IF
628# endif
629!
630! Set staring and ending observation indices.
631!
632 IF (backward) THEN
633 nstrobs(ng)=0
634 nendobs(ng)=0
635 DO i=1,ifirst
636 nstrobs(ng)=nstrobs(ng)+fourdvar(ng)%NobsSurvey(i)
637 END DO
638 nstrobs(ng)=nstrobs(ng)+1
639 ELSE
640 IF (ifirst.eq.1) THEN
641 nstrobs(ng)=0
642 nendobs(ng)=0
643 ELSE
644 nstrobs(ng)=0
645 nendobs(ng)=0
646 DO i=1,ifirst-1
647 nendobs(ng)=nendobs(ng)+fourdvar(ng)%NobsSurvey(i)
648 END DO
649 END IF
650 END IF
651!
652! Initialize total observation counter in structure array. Notice
653! that the zero index carries the total summation.
654!
655 fourdvar(ng)%ObsCount(0)=0
656 fourdvar(ng)%ObsReject(0)=0
657
658# ifdef WEAK_CONSTRAINT
659!
660! Initialize model values at observation locations.
661!
662 DO i=1,mobs
663 nlmodval(i)=inival
664 obsvetting(i)=inival
665# ifndef SP4DVAR
666 tlmodval(i)=inival
667# endif
668# if defined SP4DVAR && defined DISJOINTED && \
669 defined concurrent_kernel
670 tlmodval(i)=inival
671# endif
672 END DO
673# if defined SP4DVAR && defined DISJOINTED && \
674 !defined CONCURRENT_KERNEL
675 IF (model.ne.iadm) THEN
676 DO i=1,mobs
677 tlmodval(i)=inival
678 END DO
679 END IF
680# endif
681# endif
682!
683 10 FORMAT (/,' OBS_INITIAL_PIO - unable to open input NetCDF', &
684 & ' file: ',a)
685 20 FORMAT (/,' OBS_INITIAL_PIO - unable to find model variable:', &
686 & 1x,a,/,20x,'in input NetCDF file: ',a)
687 30 FORMAT (/,' OBS_INITIAL_PIO - No are observations available', &
688 & ' for time window (days): ',/,12x,f12.4,' - ',f12.4,/)
689!
690 RETURN
type(var_desc_t), dimension(:), pointer var_desc
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)

References mod_scalars::day2sec, mod_scalars::dt, mod_scalars::exit_flag, strings_mod::founderror(), mod_fourdvar::fourdvar, mod_fourdvar::haveobsmeta, mod_param::iadm, mod_ncparam::idnobs, mod_ncparam::idobsd, mod_ncparam::idobst, mod_ncparam::idobsx, mod_ncparam::idobsy, mod_ncparam::idobsz, mod_ncparam::idoday, mod_ncparam::idoerr, mod_ncparam::idomet, mod_ncparam::idotyp, mod_ncparam::idoval, mod_parallel::master, mod_fourdvar::mobs, mod_fourdvar::nendobs, mod_fourdvar::nlmodval, mod_fourdvar::nobsvar, mod_scalars::noerror, mod_fourdvar::nstrobs, mod_fourdvar::nsurvey, mod_scalars::ntimes, mod_scalars::ntstart, mod_ncparam::nv, mod_iounits::obs, mod_fourdvar::obssurvey, mod_scalars::obstime, mod_fourdvar::obsvetting, mod_pio_netcdf::pio_netcdf_inq_var(), mod_pio_netcdf::pio_netcdf_open(), mod_fourdvar::processobs, mod_scalars::sec2day, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::tdays, mod_scalars::time, mod_pio_netcdf::var_desc, and mod_ncparam::vname.

Referenced by obs_initial().

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