New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
diaobs.F90 in NEMO/trunk/src/OCE/OBS – NEMO

source: NEMO/trunk/src/OCE/OBS/diaobs.F90 @ 14215

Last change on this file since 14215 was 14056, checked in by ayoung, 3 years ago

Adding branch for ticket #2567 to trunk.

  • Property svn:keywords set to Id
File size: 52.3 KB
Line 
1MODULE diaobs
2   !!======================================================================
3   !!                       ***  MODULE diaobs  ***
4   !! Observation diagnostics: Computation of the misfit between data and
5   !!                          their model equivalent
6   !!======================================================================
7   !! History :  1.0  !  2006-03  (K. Mogensen) Original code
8   !!             -   !  2006-05  (K. Mogensen, A. Weaver) Reformatted
9   !!             -   !  2006-10  (A. Weaver) Cleaning and add controls
10   !!             -   !  2007-03  (K. Mogensen) General handling of profiles
11   !!             -   !  2007-04  (G. Smith) Generalized surface operators
12   !!            2.0  !  2008-10  (M. Valdivieso) obs operator for velocity profiles
13   !!            3.4  !  2014-08  (J. While) observation operator for profiles in all vertical coordinates
14   !!             -   !                      Incorporated SST bias correction 
15   !!            3.6  !  2015-02  (M. Martin) Simplification of namelist and code
16   !!             -   !  2015-08  (M. Martin) Combined surface/profile routines.
17   !!            4.0  !  2017-11  (G. Madec) style only
18   !!----------------------------------------------------------------------
19
20   !!----------------------------------------------------------------------
21   !!   dia_obs_init  : Reading and prepare observations
22   !!   dia_obs       : Compute model equivalent to observations
23   !!   dia_obs_wri   : Write observational diagnostics
24   !!   calc_date     : Compute the date of timestep in YYYYMMDD.HHMMSS format
25   !!   ini_date      : Compute the initial date YYYYMMDD.HHMMSS
26   !!   fin_date      : Compute the final date YYYYMMDD.HHMMSS
27   !!----------------------------------------------------------------------
28   USE par_kind       ! Precision variables
29   USE in_out_manager ! I/O manager
30   USE par_oce        ! ocean parameter
31   USE dom_oce        ! Ocean space and time domain variables
32   USE sbc_oce        ! Sea-ice fraction
33   !
34   USE obs_read_prof  ! Reading and allocation of profile obs
35   USE obs_read_surf  ! Reading and allocation of surface obs
36   USE obs_sstbias    ! Bias correction routine for SST
37   USE obs_readmdt    ! Reading and allocation of MDT for SLA.
38   USE obs_prep       ! Preparation of obs. (grid search etc).
39   USE obs_oper       ! Observation operators
40   USE obs_write      ! Writing of observation related diagnostics
41   USE obs_grid       ! Grid searching
42   USE obs_read_altbias ! Bias treatment for altimeter
43   USE obs_profiles_def ! Profile data definitions
44   USE obs_surf_def   ! Surface data definitions
45   USE obs_types      ! Definitions for observation types
46   !
47   USE mpp_map        ! MPP mapping
48   USE lib_mpp        ! For ctl_warn/stop
49
50   IMPLICIT NONE
51   PRIVATE
52
53   PUBLIC dia_obs_init     ! Initialize and read observations
54   PUBLIC dia_obs          ! Compute model equivalent to observations
55   PUBLIC dia_obs_wri      ! Write model equivalent to observations
56   PUBLIC dia_obs_dealloc  ! Deallocate dia_obs data
57   PUBLIC calc_date        ! Compute the date of a timestep
58
59   LOGICAL, PUBLIC :: ln_diaobs            !: Logical switch for the obs operator
60   LOGICAL         :: ln_sstnight          !  Logical switch for night mean SST obs
61   LOGICAL         :: ln_default_fp_indegs !  T=> Default obs footprint size specified in degrees, F=> in metres
62   LOGICAL         :: ln_sla_fp_indegs     !  T=> SLA obs footprint size specified in degrees, F=> in metres
63   LOGICAL         :: ln_sst_fp_indegs     !  T=> SST obs footprint size specified in degrees, F=> in metres
64   LOGICAL         :: ln_sss_fp_indegs     !  T=> SSS obs footprint size specified in degrees, F=> in metres
65   LOGICAL         :: ln_sic_fp_indegs     !  T=> sea-ice obs footprint size specified in degrees, F=> in metres
66
67   REAL(wp) ::   rn_default_avglamscl      ! E/W diameter of SLA observation footprint (metres)
68   REAL(wp) ::   rn_default_avgphiscl      ! N/S diameter of SLA observation footprint (metre
69   REAL(wp) ::   rn_sla_avglamscl          ! E/W diameter of SLA observation footprint (metres)
70   REAL(wp) ::   rn_sla_avgphiscl          ! N/S diameter of SLA observation footprint (metres)
71   REAL(wp) ::   rn_sst_avglamscl          ! E/W diameter of SST observation footprint (metres)
72   REAL(wp) ::   rn_sst_avgphiscl          ! N/S diameter of SST observation footprint (metres)
73   REAL(wp) ::   rn_sss_avglamscl          ! E/W diameter of SSS observation footprint (metres)
74   REAL(wp) ::   rn_sss_avgphiscl          ! N/S diameter of SSS observation footprint (metres)
75   REAL(wp) ::   rn_sic_avglamscl          ! E/W diameter of sea-ice observation footprint (metres)
76   REAL(wp) ::   rn_sic_avgphiscl          ! N/S diameter of sea-ice observation footprint (metres)
77
78   INTEGER :: nn_1dint                     ! Vertical interpolation method
79   INTEGER :: nn_2dint_default             ! Default horizontal interpolation method
80   INTEGER :: nn_2dint_sla                 ! SLA horizontal interpolation method
81   INTEGER :: nn_2dint_sst                 ! SST horizontal interpolation method
82   INTEGER :: nn_2dint_sss                 ! SSS horizontal interpolation method
83   INTEGER :: nn_2dint_sic                 ! Seaice horizontal interpolation method
84   INTEGER, DIMENSION(imaxavtypes) ::   nn_profdavtypes   ! Profile data types representing a daily average
85   INTEGER :: nproftypes     ! Number of profile obs types
86   INTEGER :: nsurftypes     ! Number of surface obs types
87   INTEGER , DIMENSION(:), ALLOCATABLE ::   nvarsprof, nvarssurf   ! Number of profile & surface variables
88   INTEGER , DIMENSION(:), ALLOCATABLE ::   nextrprof, nextrsurf   ! Number of profile & surface extra variables
89   INTEGER , DIMENSION(:), ALLOCATABLE ::   n2dintsurf             ! Interpolation option for surface variables
90   REAL(wp), DIMENSION(:), ALLOCATABLE ::   zavglamscl, zavgphiscl ! E/W & N/S diameter of averaging footprint for surface variables
91   LOGICAL , DIMENSION(:), ALLOCATABLE ::   lfpindegs              ! T=> surface obs footprint size specified in degrees, F=> in metres
92   LOGICAL , DIMENSION(:), ALLOCATABLE ::   llnightav              ! Logical for calculating night-time averages
93
94   TYPE(obs_surf), PUBLIC, POINTER, DIMENSION(:) ::   surfdata     !: Initial surface data
95   TYPE(obs_surf), PUBLIC, POINTER, DIMENSION(:) ::   surfdataqc   !: Surface data after quality control
96   TYPE(obs_prof), PUBLIC, POINTER, DIMENSION(:) ::   profdata     !: Initial profile data
97   TYPE(obs_prof), PUBLIC, POINTER, DIMENSION(:) ::   profdataqc   !: Profile data after quality control
98
99   CHARACTER(len=8), PUBLIC, DIMENSION(:), ALLOCATABLE ::   cobstypesprof, cobstypessurf   !: Profile & surface obs types
100
101   !!----------------------------------------------------------------------
102   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
103   !! $Id$
104   !! Software governed by the CeCILL license (see ./LICENSE)
105   !!----------------------------------------------------------------------
106CONTAINS
107
108   SUBROUTINE dia_obs_init( Kmm )
109      !!----------------------------------------------------------------------
110      !!                    ***  ROUTINE dia_obs_init  ***
111      !!         
112      !! ** Purpose : Initialize and read observations
113      !!
114      !! ** Method  : Read the namelist and call reading routines
115      !!
116      !! ** Action  : Read the namelist and call reading routines
117      !!
118      !!----------------------------------------------------------------------
119      INTEGER, INTENT(in)                ::   Kmm                      ! ocean time level indices
120      INTEGER, PARAMETER                 ::   jpmaxnfiles = 1000       ! Maximum number of files for each obs type
121      INTEGER, DIMENSION(:), ALLOCATABLE ::   ifilesprof, ifilessurf   ! Number of profile & surface files
122      INTEGER :: ios             ! Local integer output status for namelist read
123      INTEGER :: jtype           ! Counter for obs types
124      INTEGER :: jvar            ! Counter for variables
125      INTEGER :: jfile           ! Counter for files
126      INTEGER :: jnumsstbias     ! Number of SST bias files to read and apply
127      INTEGER :: n2dint_type     ! Local version of nn_2dint*
128      !
129      CHARACTER(len=128), DIMENSION(jpmaxnfiles) :: &
130         & cn_profbfiles, &      ! T/S profile input filenames
131         & cn_sstfbfiles, &      ! Sea surface temperature input filenames
132         & cn_sssfbfiles, &      ! Sea surface salinity input filenames
133         & cn_slafbfiles, &      ! Sea level anomaly input filenames
134         & cn_sicfbfiles, &      ! Seaice concentration input filenames
135         & cn_velfbfiles, &      ! Velocity profile input filenames
136         & cn_sstbiasfiles       ! SST bias input filenames
137      CHARACTER(LEN=128) :: &
138         & cn_altbiasfile        ! Altimeter bias input filename
139      CHARACTER(len=128), DIMENSION(:,:), ALLOCATABLE :: &
140         & clproffiles, &        ! Profile filenames
141         & clsurffiles           ! Surface filenames
142      CHARACTER(len=8), DIMENSION(:), ALLOCATABLE :: &
143         & clvars                ! Expected variable names
144         !
145      LOGICAL :: ln_t3d          ! Logical switch for temperature profiles
146      LOGICAL :: ln_s3d          ! Logical switch for salinity profiles
147      LOGICAL :: ln_sla          ! Logical switch for sea level anomalies
148      LOGICAL :: ln_sst          ! Logical switch for sea surface temperature
149      LOGICAL :: ln_sss          ! Logical switch for sea surface salinity
150      LOGICAL :: ln_sic          ! Logical switch for sea ice concentration
151      LOGICAL :: ln_vel3d        ! Logical switch for velocity (u,v) obs
152      LOGICAL :: ln_nea          ! Logical switch to remove obs near land
153      LOGICAL :: ln_altbias      ! Logical switch for altimeter bias
154      LOGICAL :: ln_sstbias      ! Logical switch for bias corection of SST
155      LOGICAL :: ln_ignmis       ! Logical switch for ignoring missing files
156      LOGICAL :: ln_s_at_t       ! Logical switch to compute model S at T obs
157      LOGICAL :: ln_bound_reject ! Logical to remove obs near boundaries in LAMs.
158      LOGICAL :: ltype_fp_indegs ! Local version of ln_*_fp_indegs
159      LOGICAL :: ltype_night     ! Local version of ln_sstnight (false for other variables)
160      LOGICAL, DIMENSION(:), ALLOCATABLE :: llvar   ! Logical for profile variable read
161      LOGICAL, DIMENSION(jpmaxnfiles) :: lmask ! Used for finding number of sstbias files
162      !
163      REAL(dp) :: rn_dobsini      ! Obs window start date YYYYMMDD.HHMMSS
164      REAL(dp) :: rn_dobsend      ! Obs window end date   YYYYMMDD.HHMMSS
165      REAL(wp) :: ztype_avglamscl ! Local version of rn_*_avglamscl
166      REAL(wp) :: ztype_avgphiscl ! Local version of rn_*_avgphiscl
167      REAL(wp), DIMENSION(:,:,:),   ALLOCATABLE :: zglam   ! Model longitudes for profile variables
168      REAL(wp), DIMENSION(:,:,:),   ALLOCATABLE :: zgphi   ! Model latitudes  for profile variables
169      REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: zmask   ! Model land/sea mask associated with variables
170      !!
171      NAMELIST/namobs/ln_diaobs, ln_t3d, ln_s3d, ln_sla,              &
172         &            ln_sst, ln_sic, ln_sss, ln_vel3d,               &
173         &            ln_altbias, ln_sstbias, ln_nea,                 &
174         &            ln_grid_global, ln_grid_search_lookup,          &
175         &            ln_ignmis, ln_s_at_t, ln_bound_reject,          &
176         &            ln_sstnight, ln_default_fp_indegs,              &
177         &            ln_sla_fp_indegs, ln_sst_fp_indegs,             &
178         &            ln_sss_fp_indegs, ln_sic_fp_indegs,             &
179         &            cn_profbfiles, cn_slafbfiles,                   &
180         &            cn_sstfbfiles, cn_sicfbfiles,                   &
181         &            cn_velfbfiles, cn_sssfbfiles,                   &
182         &            cn_sstbiasfiles, cn_altbiasfile,                &
183         &            cn_gridsearchfile, rn_gridsearchres,            &
184         &            rn_dobsini, rn_dobsend,                         &
185         &            rn_default_avglamscl, rn_default_avgphiscl,     &
186         &            rn_sla_avglamscl, rn_sla_avgphiscl,             &
187         &            rn_sst_avglamscl, rn_sst_avgphiscl,             &
188         &            rn_sss_avglamscl, rn_sss_avgphiscl,             &
189         &            rn_sic_avglamscl, rn_sic_avgphiscl,             &
190         &            nn_1dint, nn_2dint_default,                     &
191         &            nn_2dint_sla, nn_2dint_sst,                     &
192         &            nn_2dint_sss, nn_2dint_sic,                     &
193         &            nn_msshc, rn_mdtcorr, rn_mdtcutoff,             &
194         &            nn_profdavtypes
195      !-----------------------------------------------------------------------
196
197      !-----------------------------------------------------------------------
198      ! Read namelist parameters
199      !-----------------------------------------------------------------------
200      ! Some namelist arrays need initialising
201      cn_profbfiles  (:) = ''
202      cn_slafbfiles  (:) = ''
203      cn_sstfbfiles  (:) = ''
204      cn_sicfbfiles  (:) = ''
205      cn_velfbfiles  (:) = ''
206      cn_sssfbfiles  (:) = ''
207      cn_sstbiasfiles(:) = ''
208      nn_profdavtypes(:) = -1
209
210      CALL ini_date( rn_dobsini )
211      CALL fin_date( rn_dobsend )
212
213      ! Read namelist namobs : control observation diagnostics
214      READ  ( numnam_ref, namobs, IOSTAT = ios, ERR = 901)
215901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namobs in reference namelist' )
216      READ  ( numnam_cfg, namobs, IOSTAT = ios, ERR = 902 )
217902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namobs in configuration namelist' )
218      IF(lwm) WRITE ( numond, namobs )
219
220      IF( .NOT.ln_diaobs ) THEN
221         IF(lwp) WRITE(numout,*)
222         IF(lwp) WRITE(numout,*) 'dia_obs_init : NO Observation diagnostic used'
223         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
224         RETURN
225      ENDIF
226
227      IF(lwp) THEN
228         WRITE(numout,*)
229         WRITE(numout,*) 'dia_obs_init : Observation diagnostic initialization'
230         WRITE(numout,*) '~~~~~~~~~~~~'
231         WRITE(numout,*) '   Namelist namobs : set observation diagnostic parameters' 
232         WRITE(numout,*) '      Logical switch for T profile observations                ln_t3d = ', ln_t3d
233         WRITE(numout,*) '      Logical switch for S profile observations                ln_s3d = ', ln_s3d
234         WRITE(numout,*) '      Logical switch for SLA observations                      ln_sla = ', ln_sla
235         WRITE(numout,*) '      Logical switch for SST observations                      ln_sst = ', ln_sst
236         WRITE(numout,*) '      Logical switch for Sea Ice observations                  ln_sic = ', ln_sic
237         WRITE(numout,*) '      Logical switch for velocity observations               ln_vel3d = ', ln_vel3d
238         WRITE(numout,*) '      Logical switch for SSS observations                      ln_sss = ', ln_sss
239         WRITE(numout,*) '      Global distribution of observations              ln_grid_global = ', ln_grid_global
240         WRITE(numout,*) '      Logical switch for obs grid search lookup ln_grid_search_lookup = ', ln_grid_search_lookup
241         IF (ln_grid_search_lookup) &
242            WRITE(numout,*) '      Grid search lookup file header                cn_gridsearchfile = ', cn_gridsearchfile
243         WRITE(numout,*) '      Initial date in window YYYYMMDD.HHMMSS               rn_dobsini = ', rn_dobsini
244         WRITE(numout,*) '      Final date in window YYYYMMDD.HHMMSS                 rn_dobsend = ', rn_dobsend
245         WRITE(numout,*) '      Type of vertical interpolation method                  nn_1dint = ', nn_1dint
246         WRITE(numout,*) '      Default horizontal interpolation method        nn_2dint_default = ', nn_2dint_default
247         WRITE(numout,*) '      Type of horizontal interpolation method for SLA    nn_2dint_sla = ', nn_2dint_sla
248         WRITE(numout,*) '      Type of horizontal interpolation method for SST    nn_2dint_sst = ', nn_2dint_sst
249         WRITE(numout,*) '      Type of horizontal interpolation method for SSS    nn_2dint_sss = ', nn_2dint_sss       
250         WRITE(numout,*) '      Type of horizontal interpolation method for SIC    nn_2dint_sic = ', nn_2dint_sic
251         WRITE(numout,*) '      Default E/W diameter of obs footprint      rn_default_avglamscl = ', rn_default_avglamscl
252         WRITE(numout,*) '      Default N/S diameter of obs footprint      rn_default_avgphiscl = ', rn_default_avgphiscl
253         WRITE(numout,*) '      Default obs footprint in deg [T] or m [F]  ln_default_fp_indegs = ', ln_default_fp_indegs
254         WRITE(numout,*) '      SLA E/W diameter of obs footprint              rn_sla_avglamscl = ', rn_sla_avglamscl
255         WRITE(numout,*) '      SLA N/S diameter of obs footprint              rn_sla_avgphiscl = ', rn_sla_avgphiscl
256         WRITE(numout,*) '      SLA obs footprint in deg [T] or m [F]          ln_sla_fp_indegs = ', ln_sla_fp_indegs
257         WRITE(numout,*) '      SST E/W diameter of obs footprint              rn_sst_avglamscl = ', rn_sst_avglamscl
258         WRITE(numout,*) '      SST N/S diameter of obs footprint              rn_sst_avgphiscl = ', rn_sst_avgphiscl
259         WRITE(numout,*) '      SST obs footprint in deg [T] or m [F]          ln_sst_fp_indegs = ', ln_sst_fp_indegs
260         WRITE(numout,*) '      SIC E/W diameter of obs footprint              rn_sic_avglamscl = ', rn_sic_avglamscl
261         WRITE(numout,*) '      SIC N/S diameter of obs footprint              rn_sic_avgphiscl = ', rn_sic_avgphiscl
262         WRITE(numout,*) '      SIC obs footprint in deg [T] or m [F]          ln_sic_fp_indegs = ', ln_sic_fp_indegs
263         WRITE(numout,*) '      Rejection of observations near land switch               ln_nea = ', ln_nea
264         WRITE(numout,*) '      Rejection of obs near open bdys                 ln_bound_reject = ', ln_bound_reject
265         WRITE(numout,*) '      MSSH correction scheme                                 nn_msshc = ', nn_msshc
266         WRITE(numout,*) '      MDT  correction                                      rn_mdtcorr = ', rn_mdtcorr
267         WRITE(numout,*) '      MDT cutoff for computed correction                 rn_mdtcutoff = ', rn_mdtcutoff
268         WRITE(numout,*) '      Logical switch for alt bias                          ln_altbias = ', ln_altbias
269         WRITE(numout,*) '      Logical switch for sst bias                          ln_sstbias = ', ln_sstbias
270         WRITE(numout,*) '      Logical switch for ignoring missing files             ln_ignmis = ', ln_ignmis
271         WRITE(numout,*) '      Daily average types                             nn_profdavtypes = ', nn_profdavtypes
272         WRITE(numout,*) '      Logical switch for night-time SST obs               ln_sstnight = ', ln_sstnight
273      ENDIF
274      !-----------------------------------------------------------------------
275      ! Set up list of observation types to be used
276      ! and the files associated with each type
277      !-----------------------------------------------------------------------
278
279      nproftypes = COUNT( (/ln_t3d .OR. ln_s3d, ln_vel3d /) )
280      nsurftypes = COUNT( (/ln_sla, ln_sst, ln_sic, ln_sss /) )
281
282      IF( ln_sstbias ) THEN
283         lmask(:) = .FALSE. 
284         WHERE( cn_sstbiasfiles(:) /= '' )   lmask(:) = .TRUE. 
285         jnumsstbias = COUNT(lmask) 
286         lmask(:) = .FALSE. 
287      ENDIF     
288
289      IF( nproftypes == 0 .AND. nsurftypes == 0 ) THEN
290         CALL ctl_warn( 'dia_obs_init: ln_diaobs is set to true, but all obs operator logical flags',   &
291            &           ' (ln_t3d, ln_s3d, ln_sla, ln_sst, ln_sic, ln_vel3d)',                          &
292            &           ' are set to .FALSE. so turning off calls to dia_obs'  )
293         ln_diaobs = .FALSE.
294         RETURN
295      ENDIF
296
297      IF( nproftypes > 0 ) THEN
298         !
299         ALLOCATE( cobstypesprof(nproftypes)             )
300         ALLOCATE( ifilesprof   (nproftypes)             )
301         ALLOCATE( clproffiles  (nproftypes,jpmaxnfiles) )
302         !
303         jtype = 0
304         IF( ln_t3d .OR. ln_s3d ) THEN
305            jtype = jtype + 1
306            cobstypesprof(jtype) = 'prof'
307            clproffiles(jtype,:) = cn_profbfiles
308         ENDIF
309         IF( ln_vel3d ) THEN
310            jtype = jtype + 1
311            cobstypesprof(jtype) = 'vel'
312            clproffiles(jtype,:) = cn_velfbfiles
313         ENDIF
314         !
315         CALL obs_settypefiles( nproftypes, jpmaxnfiles, ifilesprof, cobstypesprof, clproffiles )
316         !
317      ENDIF
318
319      IF( nsurftypes > 0 ) THEN
320         !
321         ALLOCATE( cobstypessurf(nsurftypes)             )
322         ALLOCATE( ifilessurf   (nsurftypes)             )
323         ALLOCATE( clsurffiles  (nsurftypes,jpmaxnfiles) )
324         ALLOCATE( n2dintsurf   (nsurftypes)             )
325         ALLOCATE( zavglamscl   (nsurftypes)             )
326         ALLOCATE( zavgphiscl   (nsurftypes)             )
327         ALLOCATE( lfpindegs    (nsurftypes)             )
328         ALLOCATE( llnightav    (nsurftypes)             )
329         !
330         jtype = 0
331         IF( ln_sla ) THEN
332            jtype = jtype + 1
333            cobstypessurf(jtype) = 'sla'
334            clsurffiles(jtype,:) = cn_slafbfiles
335         ENDIF
336         IF( ln_sst ) THEN
337            jtype = jtype + 1
338            cobstypessurf(jtype) = 'sst'
339            clsurffiles(jtype,:) = cn_sstfbfiles
340         ENDIF
341#if defined key_si3 || defined key_cice
342         IF( ln_sic ) THEN
343            jtype = jtype + 1
344            cobstypessurf(jtype) = 'sic'
345            clsurffiles(jtype,:) = cn_sicfbfiles
346         ENDIF
347#endif
348         IF( ln_sss ) THEN
349            jtype = jtype + 1
350            cobstypessurf(jtype) = 'sss'
351            clsurffiles(jtype,:) = cn_sssfbfiles
352         ENDIF
353         !
354         CALL obs_settypefiles( nsurftypes, jpmaxnfiles, ifilessurf, cobstypessurf, clsurffiles )
355
356         DO jtype = 1, nsurftypes
357
358            IF ( TRIM(cobstypessurf(jtype)) == 'sla' ) THEN
359               IF ( nn_2dint_sla == -1 ) THEN
360                  n2dint_type  = nn_2dint_default
361               ELSE
362                  n2dint_type  = nn_2dint_sla
363               ENDIF
364               ztype_avglamscl = rn_sla_avglamscl
365               ztype_avgphiscl = rn_sla_avgphiscl
366               ltype_fp_indegs = ln_sla_fp_indegs
367               ltype_night     = .FALSE.
368            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sst' ) THEN
369               IF ( nn_2dint_sst == -1 ) THEN
370                  n2dint_type  = nn_2dint_default
371               ELSE
372                  n2dint_type  = nn_2dint_sst
373               ENDIF
374               ztype_avglamscl = rn_sst_avglamscl
375               ztype_avgphiscl = rn_sst_avgphiscl
376               ltype_fp_indegs = ln_sst_fp_indegs
377               ltype_night     = ln_sstnight
378            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sic' ) THEN
379               IF ( nn_2dint_sic == -1 ) THEN
380                  n2dint_type  = nn_2dint_default
381               ELSE
382                  n2dint_type  = nn_2dint_sic
383               ENDIF
384               ztype_avglamscl = rn_sic_avglamscl
385               ztype_avgphiscl = rn_sic_avgphiscl
386               ltype_fp_indegs = ln_sic_fp_indegs
387               ltype_night     = .FALSE.
388            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sss' ) THEN
389               IF ( nn_2dint_sss == -1 ) THEN
390                  n2dint_type  = nn_2dint_default
391               ELSE
392                  n2dint_type  = nn_2dint_sss
393               ENDIF
394               ztype_avglamscl = rn_sss_avglamscl
395               ztype_avgphiscl = rn_sss_avgphiscl
396               ltype_fp_indegs = ln_sss_fp_indegs
397               ltype_night     = .FALSE.
398            ELSE
399               n2dint_type     = nn_2dint_default
400               ztype_avglamscl = rn_default_avglamscl
401               ztype_avgphiscl = rn_default_avgphiscl
402               ltype_fp_indegs = ln_default_fp_indegs
403               ltype_night     = .FALSE.
404            ENDIF
405           
406            CALL obs_setinterpopts( nsurftypes, jtype, TRIM(cobstypessurf(jtype)), &
407               &                    nn_2dint_default, n2dint_type,                 &
408               &                    ztype_avglamscl, ztype_avgphiscl,              &
409               &                    ltype_fp_indegs, ltype_night,                  &
410               &                    n2dintsurf, zavglamscl, zavgphiscl,            &
411               &                    lfpindegs, llnightav )
412
413         END DO
414         !
415      ENDIF
416
417
418      !-----------------------------------------------------------------------
419      ! Obs operator parameter checking and initialisations
420      !-----------------------------------------------------------------------
421      !
422      IF( ln_vel3d  .AND.  .NOT.ln_grid_global ) THEN
423         CALL ctl_stop( 'Velocity data only works with ln_grid_global=.true.' )
424         RETURN
425      ENDIF
426      !
427      IF( ln_grid_global ) THEN
428         CALL ctl_warn( 'dia_obs_init: ln_grid_global=T may cause memory issues when used with a large number of processors' )
429      ENDIF
430      !
431      IF( nn_1dint < 0  .OR.  nn_1dint > 1 ) THEN
432         CALL ctl_stop('dia_obs_init: Choice of vertical (1D) interpolation method is not available')
433      ENDIF
434      !
435      IF( nn_2dint_default < 0  .OR.  nn_2dint_default > 6  ) THEN
436         CALL ctl_stop('dia_obs_init: Choice of default horizontal (2D) interpolation method is not available')
437      ENDIF
438      !
439      CALL obs_typ_init
440      IF( ln_grid_global )   CALL mppmap_init
441      !
442      CALL obs_grid_setup( )
443
444      !-----------------------------------------------------------------------
445      ! Depending on switches read the various observation types
446      !-----------------------------------------------------------------------
447      !
448      IF( nproftypes > 0 ) THEN
449         !
450         ALLOCATE( profdata  (nproftypes) , nvarsprof (nproftypes) )
451         ALLOCATE( profdataqc(nproftypes) , nextrprof (nproftypes) )
452         !
453         DO jtype = 1, nproftypes
454            !
455            IF ( TRIM(cobstypesprof(jtype)) == 'prof' ) THEN
456               nvarsprof(jtype) = 2
457               nextrprof(jtype) = 1             
458               ALLOCATE( llvar (nvarsprof(jtype)) )
459               ALLOCATE( clvars(nvarsprof(jtype)) )
460               ALLOCATE( zglam(jpi, jpj,      nvarsprof(jtype)) )
461               ALLOCATE( zgphi(jpi, jpj,      nvarsprof(jtype)) )
462               ALLOCATE( zmask(jpi, jpj, jpk, nvarsprof(jtype)) )
463               llvar(1)       = ln_t3d
464               llvar(2)       = ln_s3d
465               clvars(1)      = 'POTM'
466               clvars(2)      = 'PSAL'
467               zglam(:,:,1)   = glamt(:,:)
468               zglam(:,:,2)   = glamt(:,:)
469               zgphi(:,:,1)   = gphit(:,:)
470               zgphi(:,:,2)   = gphit(:,:)
471               zmask(:,:,:,1) = tmask(:,:,:)
472               zmask(:,:,:,2) = tmask(:,:,:)
473            ELSE IF ( TRIM(cobstypesprof(jtype)) == 'vel' )  THEN
474               nvarsprof(jtype) = 2
475               nextrprof(jtype) = 2
476               ALLOCATE( llvar (nvarsprof(jtype)) )
477               ALLOCATE( clvars(nvarsprof(jtype)) )
478               ALLOCATE( zglam(jpi, jpj,      nvarsprof(jtype)) )
479               ALLOCATE( zgphi(jpi, jpj,      nvarsprof(jtype)) )
480               ALLOCATE( zmask(jpi, jpj, jpk, nvarsprof(jtype)) )
481               llvar(1)       = ln_vel3d
482               llvar(2)       = ln_vel3d
483               clvars(1)      = 'UVEL'
484               clvars(2)      = 'VVEL'
485               zglam(:,:,1)   = glamu(:,:)
486               zglam(:,:,2)   = glamv(:,:)
487               zgphi(:,:,1)   = gphiu(:,:)
488               zgphi(:,:,2)   = gphiv(:,:)
489               zmask(:,:,:,1) = umask(:,:,:)
490               zmask(:,:,:,2) = vmask(:,:,:)
491            ELSE
492               nvarsprof(jtype) = 1
493               nextrprof(jtype) = 0
494               ALLOCATE( llvar (nvarsprof(jtype)) )
495               ALLOCATE( clvars(nvarsprof(jtype)) )
496               ALLOCATE( zglam(jpi, jpj,      nvarsprof(jtype)) )
497               ALLOCATE( zgphi(jpi, jpj,      nvarsprof(jtype)) )
498               ALLOCATE( zmask(jpi, jpj, jpk, nvarsprof(jtype)) )
499               llvar(1)       = .TRUE.
500               zglam(:,:,1)   = glamt(:,:)
501               zgphi(:,:,1)   = gphit(:,:)
502               zmask(:,:,:,1) = tmask(:,:,:)
503            ENDIF
504            !
505            ! Read in profile or profile obs types
506            CALL obs_rea_prof( profdata(jtype), ifilesprof(jtype),       &
507               &               clproffiles(jtype,1:ifilesprof(jtype)), &
508               &               nvarsprof(jtype), nextrprof(jtype), nitend-nit000+2, &
509               &               rn_dobsini, rn_dobsend, llvar, &
510               &               ln_ignmis, ln_s_at_t, .FALSE., clvars, &
511               &               kdailyavtypes = nn_profdavtypes )
512               !
513            DO jvar = 1, nvarsprof(jtype)
514               CALL obs_prof_staend( profdata(jtype), jvar )
515            END DO
516            !
517            CALL obs_pre_prof( profdata(jtype), profdataqc(jtype), &
518               &               llvar, &
519               &               jpi, jpj, jpk, &
520               &               zmask, zglam, zgphi, &
521               &               ln_nea, ln_bound_reject, Kmm, &
522               &               kdailyavtypes = nn_profdavtypes )
523            !
524            DEALLOCATE( llvar, clvars, zglam, zgphi, zmask )
525            !
526         END DO
527         !
528         DEALLOCATE( ifilesprof, clproffiles )
529         !
530      ENDIF
531      !
532      IF( nsurftypes > 0 ) THEN
533         !
534         ALLOCATE( surfdata  (nsurftypes) , nvarssurf(nsurftypes) )
535         ALLOCATE( surfdataqc(nsurftypes) , nextrsurf(nsurftypes) )
536         !
537         DO jtype = 1, nsurftypes
538            !
539            nvarssurf(jtype) = 1
540            nextrsurf(jtype) = 0
541            llnightav(jtype) = .FALSE.
542            IF( TRIM(cobstypessurf(jtype)) == 'sla' )   nextrsurf(jtype) = 2
543            IF( TRIM(cobstypessurf(jtype)) == 'sst' )   llnightav(jtype) = ln_sstnight
544            !
545            ALLOCATE( clvars( nvarssurf(jtype) ) )
546            IF ( TRIM(cobstypessurf(jtype)) == 'sla' ) THEN
547               clvars(1) = 'SLA'
548            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sst' ) THEN
549               clvars(1) = 'SST'
550            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sic' ) THEN
551               clvars(1) = 'ICECONC'
552            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sss' ) THEN
553               clvars(1) = 'SSS'
554            ENDIF
555            !
556            ! Read in surface obs types
557            CALL obs_rea_surf( surfdata(jtype), ifilessurf(jtype), &
558               &               clsurffiles(jtype,1:ifilessurf(jtype)), &
559               &               nvarssurf(jtype), nextrsurf(jtype), nitend-nit000+2, &
560               &               rn_dobsini, rn_dobsend, ln_ignmis, .FALSE., llnightav(jtype), &
561               &               clvars )
562               !
563            CALL obs_pre_surf( surfdata(jtype), surfdataqc(jtype), ln_nea, ln_bound_reject )
564            !
565            IF( TRIM(cobstypessurf(jtype)) == 'sla' ) THEN
566               CALL obs_rea_mdt( surfdataqc(jtype), n2dintsurf(jtype), Kmm )
567               IF( ln_altbias )   &
568                  & CALL obs_rea_altbias ( surfdataqc(jtype), n2dintsurf(jtype), cn_altbiasfile )
569            ENDIF
570            !
571            IF( TRIM(cobstypessurf(jtype)) == 'sst' .AND. ln_sstbias ) THEN
572               jnumsstbias = 0
573               DO jfile = 1, jpmaxnfiles
574                  IF( TRIM(cn_sstbiasfiles(jfile)) /= '' )   jnumsstbias = jnumsstbias + 1
575               END DO
576               IF( jnumsstbias == 0 )   CALL ctl_stop("ln_sstbias set but no bias files to read in")   
577               !
578               CALL obs_app_sstbias( surfdataqc(jtype), n2dintsurf(jtype)             ,   & 
579                  &                  jnumsstbias      , cn_sstbiasfiles(1:jnumsstbias) ) 
580            ENDIF
581            !
582            DEALLOCATE( clvars )
583         END DO
584         !
585         DEALLOCATE( ifilessurf, clsurffiles )
586         !
587      ENDIF
588      !
589   END SUBROUTINE dia_obs_init
590
591
592   SUBROUTINE dia_obs( kstp, Kmm )
593      !!----------------------------------------------------------------------
594      !!                    ***  ROUTINE dia_obs  ***
595      !!         
596      !! ** Purpose : Call the observation operators on each time step
597      !!
598      !! ** Method  : Call the observation operators on each time step to
599      !!              compute the model equivalent of the following data:
600      !!               - Profile data, currently T/S or U/V
601      !!               - Surface data, currently SST, SLA or sea-ice concentration.
602      !!
603      !! ** Action  :
604      !!----------------------------------------------------------------------
605      USE dom_oce, ONLY : gdept, gdept_1d     ! Ocean space domain variables (Kmm time-level only)
606      USE phycst , ONLY : rday                ! Physical constants
607      USE oce    , ONLY : ts, uu, vv, ssh     ! Ocean dynamics and tracers variables (Kmm time-level only)
608      USE phycst , ONLY : rday                ! Physical constants
609#if defined  key_si3
610      USE ice    , ONLY : at_i                ! SI3 Ice model variables
611#endif
612#if defined key_cice
613      USE sbc_oce, ONLY : fr_i     ! ice fraction
614#endif
615
616      IMPLICIT NONE
617
618      !! * Arguments
619      INTEGER, INTENT(IN) :: kstp  ! Current timestep
620      INTEGER, INTENT(in) :: Kmm   ! ocean time level indices
621      !! * Local declarations
622      INTEGER :: idaystp           ! Number of timesteps per day
623      INTEGER :: jtype             ! Data loop variable
624      INTEGER :: jvar              ! Variable number
625      INTEGER :: ji, jj            ! Loop counters
626      REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: &
627         & zprofvar                ! Model values for variables in a prof ob
628      REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: &
629         & zprofmask               ! Mask associated with zprofvar
630      REAL(wp), DIMENSION(jpi,jpj) :: &
631         & zsurfvar, &             ! Model values equivalent to surface ob.
632         & zsurfmask               ! Mask associated with surface variable
633      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: &
634         & zglam,    &             ! Model longitudes for prof variables
635         & zgphi                   ! Model latitudes for prof variables
636
637      !-----------------------------------------------------------------------
638
639      IF(lwp) THEN
640         WRITE(numout,*)
641         WRITE(numout,*) 'dia_obs : Call the observation operators', kstp
642         WRITE(numout,*) '~~~~~~~'
643      ENDIF
644
645      idaystp = NINT( rday / rn_Dt )
646
647      !-----------------------------------------------------------------------
648      ! Call the profile and surface observation operators
649      !-----------------------------------------------------------------------
650
651      IF ( nproftypes > 0 ) THEN
652
653         DO jtype = 1, nproftypes
654
655            ! Allocate local work arrays
656            ALLOCATE( zprofvar (jpi, jpj, jpk, profdataqc(jtype)%nvar) )
657            ALLOCATE( zprofmask(jpi, jpj, jpk, profdataqc(jtype)%nvar) )
658            ALLOCATE( zglam    (jpi, jpj,      profdataqc(jtype)%nvar) )
659            ALLOCATE( zgphi    (jpi, jpj,      profdataqc(jtype)%nvar) ) 
660                             
661            ! Defaults which might change
662            DO jvar = 1, profdataqc(jtype)%nvar
663               zprofmask(:,:,:,jvar) = tmask(:,:,:)
664               zglam(:,:,jvar)       = glamt(:,:)
665               zgphi(:,:,jvar)       = gphit(:,:)
666            END DO
667
668            SELECT CASE ( TRIM(cobstypesprof(jtype)) )
669            CASE('prof')
670               zprofvar(:,:,:,1) = ts(:,:,:,jp_tem,Kmm)
671               zprofvar(:,:,:,2) = ts(:,:,:,jp_sal,Kmm)
672            CASE('vel')
673               zprofvar(:,:,:,1) = uu(:,:,:,Kmm)
674               zprofvar(:,:,:,2) = vv(:,:,:,Kmm)
675               zprofmask(:,:,:,1) = umask(:,:,:)
676               zprofmask(:,:,:,2) = vmask(:,:,:)
677               zglam(:,:,1) = glamu(:,:)
678               zglam(:,:,2) = glamv(:,:)
679               zgphi(:,:,1) = gphiu(:,:)
680               zgphi(:,:,2) = gphiv(:,:)
681            CASE DEFAULT
682               CALL ctl_stop( 'Unknown profile observation type '//TRIM(cobstypesprof(jtype))//' in dia_obs' )
683            END SELECT
684
685            DO jvar = 1, profdataqc(jtype)%nvar
686               CALL obs_prof_opt( profdataqc(jtype), kstp, jpi, jpj, jpk,  &
687                  &               nit000, idaystp, jvar,                   &
688                  &               zprofvar(:,:,:,jvar),                    &
689                  &               gdept(:,:,:,Kmm), gdepw(:,:,:,Kmm),      & 
690                  &               zprofmask(:,:,:,jvar),                   &
691                  &               zglam(:,:,jvar), zgphi(:,:,jvar),        &
692                  &               nn_1dint, nn_2dint_default,              &
693                  &               kdailyavtypes = nn_profdavtypes )
694            END DO
695           
696            DEALLOCATE( zprofvar, zprofmask, zglam, zgphi )
697
698         END DO
699
700      ENDIF
701
702      IF ( nsurftypes > 0 ) THEN
703
704         DO jtype = 1, nsurftypes
705
706            !Defaults which might be changed
707            zsurfmask(:,:) = tmask(:,:,1)
708
709            SELECT CASE ( TRIM(cobstypessurf(jtype)) )
710            CASE('sst')
711               zsurfvar(:,:) = ts(:,:,1,jp_tem,Kmm)
712            CASE('sla')
713               zsurfvar(:,:) = ssh(:,:,Kmm)
714            CASE('sss')
715               zsurfvar(:,:) = ts(:,:,1,jp_sal,Kmm)
716            CASE('sic')
717               IF ( kstp == 0 ) THEN
718                  IF ( lwp .AND. surfdataqc(jtype)%nsstpmpp(1) > 0 ) THEN
719                     CALL ctl_warn( 'Sea-ice not initialised on zeroth '// &
720                        &           'time-step but some obs are valid then.' )
721                     WRITE(numout,*)surfdataqc(jtype)%nsstpmpp(1), &
722                        &           ' sea-ice obs will be missed'
723                  ENDIF
724                  surfdataqc(jtype)%nsurfup = surfdataqc(jtype)%nsurfup + &
725                     &                        surfdataqc(jtype)%nsstp(1)
726                  CYCLE
727               ELSE
728#if defined key_cice || defined key_si3
729                  zsurfvar(:,:) = fr_i(:,:)
730#else
731                  CALL ctl_stop( ' Trying to run sea-ice observation operator', &
732                     &           ' but no sea-ice model appears to have been defined' )
733#endif
734               ENDIF
735
736            END SELECT
737
738            CALL obs_surf_opt( surfdataqc(jtype), kstp, jpi, jpj,       &
739               &               nit000, idaystp, zsurfvar, zsurfmask,    &
740               &               n2dintsurf(jtype), llnightav(jtype),     &
741               &               zavglamscl(jtype), zavgphiscl(jtype),     &
742               &               lfpindegs(jtype) )
743
744         END DO
745
746      ENDIF
747
748   END SUBROUTINE dia_obs
749
750   SUBROUTINE dia_obs_wri
751      !!----------------------------------------------------------------------
752      !!                    ***  ROUTINE dia_obs_wri  ***
753      !!         
754      !! ** Purpose : Call observation diagnostic output routines
755      !!
756      !! ** Method  : Call observation diagnostic output routines
757      !!
758      !! ** Action  :
759      !!
760      !! History :
761      !!        !  06-03  (K. Mogensen) Original code
762      !!        !  06-05  (K. Mogensen) Reformatted
763      !!        !  06-10  (A. Weaver) Cleaning
764      !!        !  07-03  (K. Mogensen) General handling of profiles
765      !!        !  08-09  (M. Valdivieso) Velocity component (U,V) profiles
766      !!        !  15-08  (M. Martin) Combined writing for prof and surf types
767      !!----------------------------------------------------------------------
768      !! * Modules used
769      USE obs_rot_vel          ! Rotation of velocities
770
771      IMPLICIT NONE
772
773      !! * Local declarations
774      INTEGER :: jtype                    ! Data set loop variable
775      INTEGER :: jo, jvar, jk
776      REAL(wp), DIMENSION(:), ALLOCATABLE :: &
777         & zu, &
778         & zv
779
780      !-----------------------------------------------------------------------
781      ! Depending on switches call various observation output routines
782      !-----------------------------------------------------------------------
783
784      IF ( nproftypes > 0 ) THEN
785
786         DO jtype = 1, nproftypes
787
788            IF ( TRIM(cobstypesprof(jtype)) == 'vel' ) THEN
789
790               ! For velocity data, rotate the model velocities to N/S, E/W
791               ! using the compressed data structure.
792               ALLOCATE( &
793                  & zu(profdataqc(jtype)%nvprot(1)), &
794                  & zv(profdataqc(jtype)%nvprot(2))  &
795                  & )
796
797               CALL obs_rotvel( profdataqc(jtype), nn_2dint_default, zu, zv )
798
799               DO jo = 1, profdataqc(jtype)%nprof
800                  DO jvar = 1, 2
801                     DO jk = profdataqc(jtype)%npvsta(jo,jvar), profdataqc(jtype)%npvend(jo,jvar)
802
803                        IF ( jvar == 1 ) THEN
804                           profdataqc(jtype)%var(jvar)%vmod(jk) = zu(jk)
805                        ELSE
806                           profdataqc(jtype)%var(jvar)%vmod(jk) = zv(jk)
807                        ENDIF
808
809                     END DO
810                  END DO
811               END DO
812
813               DEALLOCATE( zu )
814               DEALLOCATE( zv )
815
816            END IF
817
818            CALL obs_prof_decompress( profdataqc(jtype), &
819               &                      profdata(jtype), .TRUE., numout )
820
821            CALL obs_wri_prof( profdata(jtype) )
822
823         END DO
824
825      ENDIF
826
827      IF ( nsurftypes > 0 ) THEN
828
829         DO jtype = 1, nsurftypes
830
831            CALL obs_surf_decompress( surfdataqc(jtype), &
832               &                      surfdata(jtype), .TRUE., numout )
833
834            CALL obs_wri_surf( surfdata(jtype) )
835
836         END DO
837
838      ENDIF
839
840   END SUBROUTINE dia_obs_wri
841
842   SUBROUTINE dia_obs_dealloc
843      IMPLICIT NONE
844      !!----------------------------------------------------------------------
845      !!                    *** ROUTINE dia_obs_dealloc ***
846      !!
847      !!  ** Purpose : To deallocate data to enable the obs_oper online loop.
848      !!               Specifically: dia_obs_init --> dia_obs --> dia_obs_wri
849      !!
850      !!  ** Method : Clean up various arrays left behind by the obs_oper.
851      !!
852      !!  ** Action :
853      !!
854      !!----------------------------------------------------------------------
855      ! obs_grid deallocation
856      CALL obs_grid_deallocate
857
858      ! diaobs deallocation
859      IF ( nproftypes > 0 ) &
860         &   DEALLOCATE( cobstypesprof, profdata, profdataqc, nvarsprof, nextrprof )
861
862      IF ( nsurftypes > 0 ) &
863         &   DEALLOCATE( cobstypessurf, surfdata, surfdataqc, nvarssurf, nextrsurf, &
864         &               n2dintsurf, zavglamscl, zavgphiscl, lfpindegs, llnightav )
865
866   END SUBROUTINE dia_obs_dealloc
867
868   SUBROUTINE calc_date( kstp, ddobs )
869      !!----------------------------------------------------------------------
870      !!                    ***  ROUTINE calc_date  ***
871      !!         
872      !! ** Purpose : Get date in double precision YYYYMMDD.HHMMSS format
873      !!
874      !! ** Method  : Get date in double precision YYYYMMDD.HHMMSS format
875      !!
876      !! ** Action  : Get date in double precision YYYYMMDD.HHMMSS format
877      !!
878      !! ** Action  : Get initial date in double precision YYYYMMDD.HHMMSS format
879      !!
880      !! History :
881      !!        !  06-03  (K. Mogensen)  Original code
882      !!        !  06-05  (K. Mogensen)  Reformatted
883      !!        !  06-10  (A. Weaver) Cleaning
884      !!        !  06-10  (G. Smith) Calculates initial date the same as method for final date
885      !!        !  10-05  (D. Lea) Update to month length calculation for NEMO vn3.2
886      !!        !  2014-09  (D. Lea) New generic routine now deals with arbitrary initial time of day
887      !!----------------------------------------------------------------------
888      USE phycst, ONLY : &            ! Physical constants
889         & rday
890      USE dom_oce, ONLY : &           ! Ocean space and time domain variables
891         & rn_Dt
892
893      IMPLICIT NONE
894
895      !! * Arguments
896      REAL(KIND=dp), INTENT(OUT) :: ddobs                        ! Date in YYYYMMDD.HHMMSS
897      INTEGER :: kstp
898
899      !! * Local declarations
900      INTEGER :: iyea        ! date - (year, month, day, hour, minute)
901      INTEGER :: imon
902      INTEGER :: iday
903      INTEGER :: ihou
904      INTEGER :: imin
905      INTEGER :: imday       ! Number of days in month.
906      REAL(wp) :: zdayfrc    ! Fraction of day
907
908      INTEGER, DIMENSION(12) ::   imonth_len    !: length in days of the months of the current year
909
910      !!----------------------------------------------------------------------
911      !! Initial date initialization (year, month, day, hour, minute)
912      !!----------------------------------------------------------------------
913      iyea =   ndate0 / 10000
914      imon = ( ndate0 - iyea * 10000 ) / 100
915      iday =   ndate0 - iyea * 10000 - imon * 100
916      ihou =   nn_time0 / 100
917      imin = ( nn_time0 - ihou * 100 ) 
918
919      !!----------------------------------------------------------------------
920      !! Compute number of days + number of hours + min since initial time
921      !!----------------------------------------------------------------------
922      zdayfrc = kstp * rn_Dt / rday
923      zdayfrc = zdayfrc - aint(zdayfrc)
924      imin = imin + int( zdayfrc * 24 * 60 ) 
925      DO WHILE (imin >= 60) 
926        imin=imin-60
927        ihou=ihou+1
928      END DO
929      DO WHILE (ihou >= 24)
930        ihou=ihou-24
931        iday=iday+1
932      END DO
933      iday = iday + kstp * rn_Dt / rday 
934
935      !-----------------------------------------------------------------------
936      ! Convert number of days (iday) into a real date
937      !----------------------------------------------------------------------
938
939      CALL calc_month_len( iyea, imonth_len )
940
941      DO WHILE ( iday > imonth_len(imon) )
942         iday = iday - imonth_len(imon)
943         imon = imon + 1 
944         IF ( imon > 12 ) THEN
945            imon = 1
946            iyea = iyea + 1
947            CALL calc_month_len( iyea, imonth_len )  ! update month lengths
948         ENDIF
949      END DO
950
951      !----------------------------------------------------------------------
952      ! Convert it into YYYYMMDD.HHMMSS format.
953      !----------------------------------------------------------------------
954      ddobs = iyea * 10000_dp + imon * 100_dp + &
955         &    iday + ihou * 0.01_dp + imin * 0.0001_dp
956
957   END SUBROUTINE calc_date
958
959   SUBROUTINE ini_date( ddobsini )
960      !!----------------------------------------------------------------------
961      !!                    ***  ROUTINE ini_date  ***
962      !!         
963      !! ** Purpose : Get initial date in double precision YYYYMMDD.HHMMSS format
964      !!
965      !! ** Method  :
966      !!
967      !! ** Action  :
968      !!
969      !! History :
970      !!        !  06-03  (K. Mogensen)  Original code
971      !!        !  06-05  (K. Mogensen)  Reformatted
972      !!        !  06-10  (A. Weaver) Cleaning
973      !!        !  10-05  (D. Lea) Update to month length calculation for NEMO vn3.2
974      !!        !  2014-09  (D. Lea) Change to call generic routine calc_date
975      !!----------------------------------------------------------------------
976
977      IMPLICIT NONE
978
979      !! * Arguments
980      REAL(KIND=dp), INTENT(OUT) :: ddobsini                   ! Initial date in YYYYMMDD.HHMMSS
981
982      CALL calc_date( nit000 - 1, ddobsini )
983
984   END SUBROUTINE ini_date
985
986   SUBROUTINE fin_date( ddobsfin )
987      !!----------------------------------------------------------------------
988      !!                    ***  ROUTINE fin_date  ***
989      !!         
990      !! ** Purpose : Get final date in double precision YYYYMMDD.HHMMSS format
991      !!
992      !! ** Method  :
993      !!
994      !! ** Action  :
995      !!
996      !! History :
997      !!        !  06-03  (K. Mogensen)  Original code
998      !!        !  06-05  (K. Mogensen)  Reformatted
999      !!        !  06-10  (A. Weaver) Cleaning
1000      !!        !  10-05  (D. Lea) Update to month length calculation for NEMO vn3.2
1001      !!        !  2014-09  (D. Lea) Change to call generic routine calc_date
1002      !!----------------------------------------------------------------------
1003
1004      IMPLICIT NONE
1005
1006      !! * Arguments
1007      REAL(dp), INTENT(OUT) :: ddobsfin ! Final date in YYYYMMDD.HHMMSS
1008
1009      CALL calc_date( nitend, ddobsfin )
1010
1011   END SUBROUTINE fin_date
1012   
1013   SUBROUTINE obs_settypefiles( ntypes, jpmaxnfiles, ifiles, cobstypes, cfiles )
1014
1015      INTEGER, INTENT(IN) :: ntypes      ! Total number of obs types
1016      INTEGER, INTENT(IN) :: jpmaxnfiles ! Maximum number of files allowed for each type
1017      INTEGER, DIMENSION(ntypes), INTENT(OUT) :: &
1018         &                   ifiles      ! Out number of files for each type
1019      CHARACTER(len=8), DIMENSION(ntypes), INTENT(IN) :: &
1020         &                   cobstypes   ! List of obs types
1021      CHARACTER(len=128), DIMENSION(ntypes, jpmaxnfiles), INTENT(IN) :: &
1022         &                   cfiles      ! List of files for all types
1023
1024      !Local variables
1025      INTEGER :: jfile
1026      INTEGER :: jtype
1027
1028      DO jtype = 1, ntypes
1029
1030         ifiles(jtype) = 0
1031         DO jfile = 1, jpmaxnfiles
1032            IF ( trim(cfiles(jtype,jfile)) /= '' ) &
1033                      ifiles(jtype) = ifiles(jtype) + 1
1034         END DO
1035
1036         IF ( ifiles(jtype) == 0 ) THEN
1037              CALL ctl_stop( 'Logical for observation type '//TRIM(cobstypes(jtype))//   &
1038                 &           ' set to true but no files available to read' )
1039         ENDIF
1040
1041         IF(lwp) THEN   
1042            WRITE(numout,*) '             '//cobstypes(jtype)//' input observation file names:'
1043            DO jfile = 1, ifiles(jtype)
1044               WRITE(numout,*) '                '//TRIM(cfiles(jtype,jfile))
1045            END DO
1046         ENDIF
1047
1048      END DO
1049
1050   END SUBROUTINE obs_settypefiles
1051
1052   SUBROUTINE obs_setinterpopts( ntypes, jtype, ctypein,             &
1053              &                  n2dint_default, n2dint_type,        &
1054              &                  ravglamscl_type, ravgphiscl_type,   &
1055              &                  lfp_indegs_type, lavnight_type,     &
1056              &                  n2dint, ravglamscl, ravgphiscl,     &
1057              &                  lfpindegs, lavnight )
1058
1059      INTEGER, INTENT(IN)  :: ntypes             ! Total number of obs types
1060      INTEGER, INTENT(IN)  :: jtype              ! Index of the current type of obs
1061      INTEGER, INTENT(IN)  :: n2dint_default     ! Default option for interpolation type
1062      INTEGER, INTENT(IN)  :: n2dint_type        ! Option for interpolation type
1063      REAL(wp), INTENT(IN) :: &
1064         &                    ravglamscl_type, & !E/W diameter of obs footprint for this type
1065         &                    ravgphiscl_type    !N/S diameter of obs footprint for this type
1066      LOGICAL, INTENT(IN)  :: lfp_indegs_type    !T=> footprint in degrees, F=> in metres
1067      LOGICAL, INTENT(IN)  :: lavnight_type      !T=> obs represent night time average
1068      CHARACTER(len=8), INTENT(IN) :: ctypein 
1069
1070      INTEGER, DIMENSION(ntypes), INTENT(INOUT) :: &
1071         &                    n2dint 
1072      REAL(wp), DIMENSION(ntypes), INTENT(INOUT) :: &
1073         &                    ravglamscl, ravgphiscl
1074      LOGICAL, DIMENSION(ntypes), INTENT(INOUT) :: &
1075         &                    lfpindegs, lavnight
1076
1077      lavnight(jtype) = lavnight_type
1078
1079      IF ( (n2dint_type >= 0) .AND. (n2dint_type <= 6) ) THEN
1080         n2dint(jtype) = n2dint_type
1081      ELSE IF ( n2dint_type == -1 ) THEN
1082         n2dint(jtype) = n2dint_default
1083      ELSE
1084         CALL ctl_stop(' Choice of '//TRIM(ctypein)//' horizontal (2D) interpolation method', &
1085           &                    ' is not available')
1086      ENDIF
1087
1088      ! For averaging observation footprints set options for size of footprint
1089      IF ( (n2dint(jtype) > 4) .AND. (n2dint(jtype) <= 6) ) THEN
1090         IF ( ravglamscl_type > 0._wp ) THEN
1091            ravglamscl(jtype) = ravglamscl_type
1092         ELSE
1093            CALL ctl_stop( 'Incorrect value set for averaging footprint '// &
1094                           'scale (ravglamscl) for observation type '//TRIM(ctypein) )     
1095         ENDIF
1096
1097         IF ( ravgphiscl_type > 0._wp ) THEN
1098            ravgphiscl(jtype) = ravgphiscl_type
1099         ELSE
1100            CALL ctl_stop( 'Incorrect value set for averaging footprint '// &
1101                           'scale (ravgphiscl) for observation type '//TRIM(ctypein) )     
1102         ENDIF
1103
1104         lfpindegs(jtype) = lfp_indegs_type 
1105
1106      ENDIF
1107
1108      ! Write out info
1109      IF(lwp) THEN
1110         IF ( n2dint(jtype) <= 4 ) THEN
1111            WRITE(numout,*) '             '//TRIM(ctypein)// &
1112               &            ' model counterparts will be interpolated horizontally'
1113         ELSE IF ( n2dint(jtype) <= 6 ) THEN
1114            WRITE(numout,*) '             '//TRIM(ctypein)// &
1115               &            ' model counterparts will be averaged horizontally'
1116            WRITE(numout,*) '             '//'    with E/W scale: ',ravglamscl(jtype)
1117            WRITE(numout,*) '             '//'    with N/S scale: ',ravgphiscl(jtype)
1118            IF ( lfpindegs(jtype) ) THEN
1119                WRITE(numout,*) '             '//'    (in degrees)'
1120            ELSE
1121                WRITE(numout,*) '             '//'    (in metres)'
1122            ENDIF
1123         ENDIF
1124      ENDIF
1125
1126   END SUBROUTINE obs_setinterpopts
1127
1128END MODULE diaobs
Note: See TracBrowser for help on using the repository browser.