obs_read_sla.F90 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548
  1. MODULE obs_read_sla
  2. !!======================================================================
  3. !! *** MODULE obs_read_sla ***
  4. !! Observation diagnostics: Read the along track SLA data from
  5. !! the AVISO/CLS database or any SLA data
  6. !! from feedback files
  7. !!======================================================================
  8. !!----------------------------------------------------------------------
  9. !! obs_rea_sla : Driver for reading SLA data from the AVISO/CLS/feedback
  10. !!----------------------------------------------------------------------
  11. !! * Modules used
  12. USE par_kind ! Precision variables
  13. USE in_out_manager ! I/O manager
  14. USE dom_oce ! Ocean space and time domain variables
  15. USE obs_mpp ! MPP support routines for observation diagnostics
  16. USE julian ! Julian date routines
  17. USE obs_utils ! Observation operator utility functions
  18. USE obs_grid ! Grid search
  19. USE obs_sort ! Sorting observation arrays
  20. USE obs_surf_def ! Surface observation definitions
  21. USE obs_types ! Observation type definitions
  22. USE obs_sla_io ! I/O for sla files
  23. USE netcdf ! NetCDF library
  24. IMPLICIT NONE
  25. !! * Routine accessibility
  26. PRIVATE
  27. PUBLIC obs_rea_sla ! Read the SLA observations from the AVISO/SLA database
  28. !!----------------------------------------------------------------------
  29. !! NEMO/OPA 3.3 , NEMO Consortium (2010)
  30. !! $Id: obs_read_sla.F90 4990 2014-12-15 16:42:49Z timgraham $
  31. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  32. !!----------------------------------------------------------------------
  33. CONTAINS
  34. SUBROUTINE obs_rea_sla( kformat, &
  35. & sladata, knumfiles, cfilenames, &
  36. & kvars, kextr, kstp, ddobsini, ddobsend, &
  37. & ldignmis, ldmod, ldobstd )
  38. !!---------------------------------------------------------------------
  39. !!
  40. !! *** ROUTINE obs_rea_sla ***
  41. !!
  42. !! ** Purpose : Read from file the SLA data
  43. !!
  44. !! ** Method : Depending on kformat either AVISO or
  45. !! feedback data files are read
  46. !!
  47. !! ** Action :
  48. !!
  49. !!
  50. !! History :
  51. !! ! : 2009-01 (K. Mogensen) Initial version based on old versions
  52. !!----------------------------------------------------------------------
  53. !! * Modules used
  54. !! * Arguments
  55. INTEGER :: kformat ! Format of input data
  56. ! ! 0: Feedback
  57. ! ! 1: AVISO
  58. TYPE(obs_surf), INTENT(INOUT) :: sladata ! SLA data to be read
  59. INTEGER, INTENT(IN) :: knumfiles ! Number of files to read in
  60. CHARACTER(LEN=128), INTENT(IN) :: &
  61. & cfilenames(knumfiles) ! File names to read in
  62. INTEGER, INTENT(IN) :: kvars ! Number of variables in sladata
  63. INTEGER, INTENT(IN) :: kextr ! Number of extra fields for each var in sladata
  64. INTEGER, INTENT(IN) :: kstp ! Ocean time-step index
  65. LOGICAL, INTENT(IN) :: ldignmis ! Ignore missing files
  66. LOGICAL, INTENT(IN) :: ldmod ! Initialize model from input data
  67. LOGICAL, INTENT(INOUT), optional :: &
  68. & ldobstd ! Read observation standard deviation from fb. file
  69. REAL(KIND=dp), INTENT(IN) :: ddobsini ! Obs. ini time in YYYYMMDD.HHMMSS
  70. REAL(KIND=dp), INTENT(IN) :: ddobsend ! Obs. end time in YYYYMMDD.HHMMSS
  71. !! * Local declarations
  72. CHARACTER(LEN=11), PARAMETER :: cpname='obs_rea_sla'
  73. INTEGER :: ji
  74. INTEGER :: jj
  75. INTEGER :: jk
  76. INTEGER :: iflag
  77. INTEGER :: inobf
  78. INTEGER :: i_file_id
  79. INTEGER :: inowin
  80. INTEGER :: iyea
  81. INTEGER :: imon
  82. INTEGER :: iday
  83. INTEGER :: ihou
  84. INTEGER :: imin
  85. INTEGER :: isec
  86. INTEGER, DIMENSION(knumfiles) :: &
  87. & irefdate
  88. INTEGER :: iobsmpp
  89. INTEGER, DIMENSION(imaxmissions+1) :: &
  90. & ityp, &
  91. & itypmpp
  92. INTEGER, DIMENSION(:), ALLOCATABLE :: &
  93. & iobsi, &
  94. & iobsj, &
  95. & iproc, &
  96. & iindx, &
  97. & ifileidx, &
  98. & islaidx
  99. INTEGER :: itype
  100. REAL(wp), DIMENSION(:), ALLOCATABLE :: &
  101. & zphi, &
  102. & zlam
  103. real(wp), DIMENSION(:), ALLOCATABLE :: &
  104. & zdat
  105. LOGICAL :: llvalprof
  106. LOGICAL :: llobstd
  107. TYPE(obfbdata), POINTER, DIMENSION(:) :: &
  108. & inpfiles
  109. real(wp), DIMENSION(knumfiles) :: &
  110. & djulini, &
  111. & djulend
  112. INTEGER, DIMENSION(knumfiles) :: &
  113. & iobspos, &
  114. & iobcpos
  115. INTEGER :: iobs
  116. INTEGER :: iobstot
  117. INTEGER :: ios
  118. INTEGER :: ioserrcount
  119. CHARACTER(len=8) :: cl_refdate
  120. ! Local initialization
  121. iobs = 0
  122. IF ( PRESENT(ldobstd) ) THEN
  123. IF (.NOT.ldmod) THEN
  124. llobstd = .false.
  125. ELSE
  126. llobstd = ldobstd
  127. ENDIF
  128. ELSE
  129. llobstd = .FALSE.
  130. ENDIF
  131. !-----------------------------------------------------------------------
  132. ! Check that the model part is just with feedback data files
  133. !-----------------------------------------------------------------------
  134. IF ( ldmod .AND. ( kformat /= 0 ) ) THEN
  135. CALL ctl_stop( 'Model can only be read from feedback data' )
  136. RETURN
  137. ENDIF
  138. !-----------------------------------------------------------------------
  139. ! Check that the prescribed obs err is just with feedback data files
  140. !-----------------------------------------------------------------------
  141. IF ( llobstd .AND. ( kformat /= 0 ) ) THEN
  142. CALL ctl_stop( 'Observation error can only be read from feedback files' )
  143. RETURN
  144. ENDIF
  145. !-----------------------------------------------------------------------
  146. ! Count the number of files needed and allocate the obfbdata type
  147. !-----------------------------------------------------------------------
  148. inobf = knumfiles
  149. ALLOCATE( inpfiles(inobf) )
  150. sla_files : DO jj = 1, inobf
  151. !---------------------------------------------------------------------
  152. ! Prints
  153. !---------------------------------------------------------------------
  154. IF(lwp) THEN
  155. WRITE(numout,*)
  156. WRITE(numout,*) ' obs_rea_sla : Reading from file = ', &
  157. & TRIM( TRIM( cfilenames(jj) ) )
  158. WRITE(numout,*) ' ~~~~~~~~~~~'
  159. WRITE(numout,*)
  160. ENDIF
  161. !---------------------------------------------------------------------
  162. ! Initialization: Open file and get dimensions only
  163. !---------------------------------------------------------------------
  164. iflag = nf90_open( TRIM( TRIM( cfilenames(jj) ) ), nf90_nowrite, &
  165. & i_file_id )
  166. IF ( iflag /= nf90_noerr ) THEN
  167. IF ( ldignmis ) THEN
  168. inpfiles(jj)%nobs = 0
  169. CALL ctl_warn( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // &
  170. & ' not found' )
  171. ELSE
  172. CALL ctl_stop( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // &
  173. & ' not found' )
  174. ENDIF
  175. ELSE
  176. !------------------------------------------------------------------
  177. ! Close the file since it is opened in read_proffile
  178. !------------------------------------------------------------------
  179. iflag = nf90_close( i_file_id )
  180. !------------------------------------------------------------------
  181. ! Read the profile file into inpfiles
  182. !------------------------------------------------------------------
  183. IF ( kformat == 0 ) THEN
  184. CALL init_obfbdata( inpfiles(jj) )
  185. IF(lwp) THEN
  186. WRITE(numout,*)
  187. WRITE(numout,*)'Reading from feedback file :', &
  188. & TRIM( cfilenames(jj) )
  189. ENDIF
  190. CALL read_obfbdata( TRIM( cfilenames(jj) ), inpfiles(jj), &
  191. & ldgrid = .TRUE. )
  192. IF ( inpfiles(jj)%nvar < 1 ) THEN
  193. CALL ctl_stop( 'Feedback format error' )
  194. RETURN
  195. ENDIF
  196. IF ( TRIM(inpfiles(jj)%cname(1)) /= 'SLA' ) THEN
  197. CALL ctl_stop( 'Feedback format error' )
  198. RETURN
  199. ENDIF
  200. IF ( ldmod .AND. ( ( inpfiles(jj)%nadd == 0 ) .OR.&
  201. & ( inpfiles(jj)%next == 0 ) ) ) THEN
  202. CALL ctl_stop( 'Model not in input data' )
  203. RETURN
  204. ENDIF
  205. ELSEIF ( kformat == 1) THEN
  206. CALL read_avisofile( TRIM( cfilenames(jj) ), inpfiles(jj), &
  207. & numout, lwp, .TRUE. )
  208. ELSE
  209. CALL ctl_stop( 'File format unknown' )
  210. ENDIF
  211. !------------------------------------------------------------------
  212. ! Change longitude (-180,180)
  213. !------------------------------------------------------------------
  214. DO ji = 1, inpfiles(jj)%nobs
  215. IF ( inpfiles(jj)%plam(ji) < -180. ) &
  216. & inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) + 360.
  217. IF ( inpfiles(jj)%plam(ji) > 180. ) &
  218. & inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) - 360.
  219. END DO
  220. !------------------------------------------------------------------
  221. ! Calculate the date (change eventually)
  222. !------------------------------------------------------------------
  223. cl_refdate=inpfiles(jj)%cdjuldref(1:8)
  224. READ(cl_refdate,'(I8)') irefdate(jj)
  225. CALL ddatetoymdhms( ddobsini, iyea, imon, iday, ihou, imin, isec )
  226. CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulini(jj), &
  227. & krefdate = irefdate(jj) )
  228. CALL ddatetoymdhms( ddobsend, iyea, imon, iday, ihou, imin, isec )
  229. CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulend(jj), &
  230. & krefdate = irefdate(jj) )
  231. IF ( inpfiles(jj)%nobs > 0 ) THEN
  232. inpfiles(jj)%iproc = -1
  233. inpfiles(jj)%iobsi = -1
  234. inpfiles(jj)%iobsj = -1
  235. ENDIF
  236. inowin = 0
  237. DO ji = 1, inpfiles(jj)%nobs
  238. IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
  239. IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
  240. IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. &
  241. & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
  242. inowin = inowin + 1
  243. ENDIF
  244. END DO
  245. ALLOCATE( zlam(inowin) )
  246. ALLOCATE( zphi(inowin) )
  247. ALLOCATE( iobsi(inowin) )
  248. ALLOCATE( iobsj(inowin) )
  249. ALLOCATE( iproc(inowin) )
  250. inowin = 0
  251. DO ji = 1, inpfiles(jj)%nobs
  252. IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
  253. IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
  254. IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. &
  255. & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
  256. inowin = inowin + 1
  257. zlam(inowin) = inpfiles(jj)%plam(ji)
  258. zphi(inowin) = inpfiles(jj)%pphi(ji)
  259. ENDIF
  260. END DO
  261. CALL obs_grid_search( inowin, zlam, zphi, iobsi, iobsj, iproc, 'T' )
  262. inowin = 0
  263. DO ji = 1, inpfiles(jj)%nobs
  264. IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
  265. IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
  266. IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. &
  267. & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
  268. inowin = inowin + 1
  269. inpfiles(jj)%iproc(ji,1) = iproc(inowin)
  270. inpfiles(jj)%iobsi(ji,1) = iobsi(inowin)
  271. inpfiles(jj)%iobsj(ji,1) = iobsj(inowin)
  272. ENDIF
  273. END DO
  274. DEALLOCATE( zlam, zphi, iobsi, iobsj, iproc )
  275. DO ji = 1, inpfiles(jj)%nobs
  276. IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
  277. IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
  278. IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. &
  279. & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
  280. IF ( nproc == 0 ) THEN
  281. IF ( inpfiles(jj)%iproc(ji,1) > nproc ) CYCLE
  282. ELSE
  283. IF ( inpfiles(jj)%iproc(ji,1) /= nproc ) CYCLE
  284. ENDIF
  285. llvalprof = .FALSE.
  286. IF ( ( inpfiles(jj)%ivlqc(1,ji,1) == 1 ) .OR. &
  287. & ( inpfiles(jj)%ivlqc(1,ji,1) == 2 ) ) THEN
  288. iobs = iobs + 1
  289. ENDIF
  290. ENDIF
  291. END DO
  292. ENDIF
  293. END DO sla_files
  294. IF (llobstd) THEN
  295. DO jj = 1, inobf
  296. iobspos(jj) = -1
  297. iobcpos(jj) = -1
  298. DO ji = 1,inpfiles(jj)%nadd
  299. IF ( TRIM(inpfiles(jj)%caddname(ji)) == 'OSTD' ) THEN
  300. iobspos(jj)=ji
  301. ENDIF
  302. IF ( TRIM(inpfiles(jj)%caddname(ji)) == 'OCNT' ) THEN
  303. iobcpos(jj)=ji
  304. ENDIF
  305. END DO
  306. END DO
  307. llobstd = ( ( MINVAL(iobspos) > 0 ) .AND. ( MINVAL(iobcpos) > 0 ) )
  308. IF (llobstd) THEN
  309. IF (lwp) WRITE(numout,*)'SLA superobs information present.'
  310. ELSE
  311. IF (lwp) WRITE(numout,*)'SLA superobs information not present.'
  312. ENDIF
  313. ENDIF
  314. !-----------------------------------------------------------------------
  315. ! Get the time ordered indices of the input data
  316. !-----------------------------------------------------------------------
  317. !---------------------------------------------------------------------
  318. ! Loop over input data files to count total number of profiles
  319. !---------------------------------------------------------------------
  320. iobstot = 0
  321. DO jj = 1, inobf
  322. DO ji = 1, inpfiles(jj)%nobs
  323. IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
  324. IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
  325. IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. &
  326. & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
  327. iobstot = iobstot + 1
  328. ENDIF
  329. END DO
  330. END DO
  331. ALLOCATE( iindx(iobstot), ifileidx(iobstot), &
  332. & islaidx(iobstot), zdat(iobstot) )
  333. jk = 0
  334. DO jj = 1, inobf
  335. DO ji = 1, inpfiles(jj)%nobs
  336. IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
  337. IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
  338. IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. &
  339. & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
  340. jk = jk + 1
  341. ifileidx(jk) = jj
  342. islaidx(jk) = ji
  343. zdat(jk) = inpfiles(jj)%ptim(ji)
  344. ENDIF
  345. END DO
  346. END DO
  347. CALL sort_dp_indx( iobstot, &
  348. & zdat, &
  349. & iindx )
  350. CALL obs_surf_alloc( sladata, iobs, kvars, kextr, kstp, jpi, jpj )
  351. ! * Read obs/positions, QC, all variable and assign to sladata
  352. iobs = 0
  353. ityp (:) = 0
  354. itypmpp(:) = 0
  355. ioserrcount = 0
  356. DO jk = 1, iobstot
  357. jj = ifileidx(iindx(jk))
  358. ji = islaidx(iindx(jk))
  359. IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
  360. IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
  361. IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. &
  362. & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
  363. IF ( nproc == 0 ) THEN
  364. IF ( inpfiles(jj)%iproc(ji,1) > nproc ) CYCLE
  365. ELSE
  366. IF ( inpfiles(jj)%iproc(ji,1) /= nproc ) CYCLE
  367. ENDIF
  368. ! Set observation information
  369. IF ( ( inpfiles(jj)%ivlqc(1,ji,1) == 1 ) .OR. &
  370. & ( inpfiles(jj)%ivlqc(1,ji,1) == 2 ) ) THEN
  371. iobs = iobs + 1
  372. CALL jul2greg( isec, &
  373. & imin, &
  374. & ihou, &
  375. & iday, &
  376. & imon, &
  377. & iyea, &
  378. & inpfiles(jj)%ptim(ji), &
  379. & irefdate(jj) )
  380. ! SLA time coordinates
  381. sladata%nyea(iobs) = iyea
  382. sladata%nmon(iobs) = imon
  383. sladata%nday(iobs) = iday
  384. sladata%nhou(iobs) = ihou
  385. sladata%nmin(iobs) = imin
  386. ! SLA space coordinates
  387. sladata%rlam(iobs) = inpfiles(jj)%plam(ji)
  388. sladata%rphi(iobs) = inpfiles(jj)%pphi(ji)
  389. ! Coordinate search parameters
  390. sladata%mi (iobs) = inpfiles(jj)%iobsi(ji,1)
  391. sladata%mj (iobs) = inpfiles(jj)%iobsj(ji,1)
  392. ! Instrument type
  393. READ( inpfiles(jj)%cdtyp(ji), '(I4)', IOSTAT = ios, ERR = 901 ) itype
  394. 901 IF ( ios /= 0 ) THEN
  395. IF (ioserrcount == 0) CALL ctl_warn ( 'Problem converting an instrument type to integer. Setting type to zero' )
  396. ioserrcount = ioserrcount + 1
  397. itype = 0
  398. ENDIF
  399. sladata%ntyp(iobs) = itype
  400. ityp(itype+1) = ityp(itype+1) + 1
  401. ! Identifier
  402. sladata%cwmo(iobs) = inpfiles(jj)%cdwmo(ji)
  403. ! Bookkeeping data to match observations
  404. sladata%nsidx(iobs) = iobs
  405. sladata%nsfil(iobs) = iindx(jk)
  406. ! QC flags
  407. sladata%nqc(iobs) = inpfiles(jj)%ivqc(ji,1)
  408. ! Observed value
  409. sladata%robs(iobs,1) = inpfiles(jj)%pob(1,ji,1)
  410. ! Model and MDT is set to fbrmdi unless read from file
  411. IF ( ldmod ) THEN
  412. sladata%rmod(iobs,1) = inpfiles(jj)%padd(1,ji,1,1)
  413. sladata%rext(iobs,1) = inpfiles(jj)%padd(1,ji,2,1)
  414. sladata%rext(iobs,2) = inpfiles(jj)%pext(1,ji,1)
  415. IF (llobstd) THEN
  416. sladata%rext(iobs,3) = &
  417. & inpfiles(jj)%padd(1,ji,iobspos(jj),1)
  418. sladata%rext(iobs,4) = &
  419. & inpfiles(jj)%padd(1,ji,iobcpos(jj),1)
  420. ENDIF
  421. ELSE
  422. sladata%rmod(iobs,1) = fbrmdi
  423. sladata%rext(iobs,:) = fbrmdi
  424. ENDIF
  425. ENDIF
  426. ENDIF
  427. END DO
  428. !-----------------------------------------------------------------------
  429. ! Sum up over processors
  430. !-----------------------------------------------------------------------
  431. CALL obs_mpp_sum_integer( iobs, iobsmpp )
  432. CALL obs_mpp_sum_integers( ityp, itypmpp, imaxmissions + 1 )
  433. !-----------------------------------------------------------------------
  434. ! Output number of observations.
  435. !-----------------------------------------------------------------------
  436. IF (lwp) THEN
  437. WRITE(numout,*)
  438. WRITE(numout,'(1X,A)')'Altimeter satellites'
  439. WRITE(numout,'(1X,A)')'--------------------'
  440. DO jj = 1,8
  441. IF ( itypmpp(jj) > 0 ) THEN
  442. WRITE(numout,'(1X,A38,A2,I10)')calttyp(jj-1),'= ',itypmpp(jj)
  443. ENDIF
  444. END DO
  445. WRITE(numout,'(1X,A50)')'--------------------------------------------------'
  446. WRITE(numout,'(1X,A40,I10)')'Total = ',iobsmpp
  447. WRITE(numout,*)
  448. ENDIF
  449. !-----------------------------------------------------------------------
  450. ! Deallocate temporary data
  451. !-----------------------------------------------------------------------
  452. DEALLOCATE( ifileidx, islaidx, zdat )
  453. !-----------------------------------------------------------------------
  454. ! Deallocate input data
  455. !-----------------------------------------------------------------------
  456. DO jj = 1, inobf
  457. CALL dealloc_obfbdata( inpfiles(jj) )
  458. END DO
  459. DEALLOCATE( inpfiles )
  460. !-----------------------------------------------------------------------
  461. ! Reset ldobstd if the data is present
  462. !-----------------------------------------------------------------------
  463. IF ( PRESENT(ldobstd) ) THEN
  464. ldobstd = llobstd
  465. ENDIF
  466. END SUBROUTINE obs_rea_sla
  467. END MODULE obs_read_sla