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/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/OBS – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/OBS/diaobs.F90 @ 14219

Last change on this file since 14219 was 14219, checked in by mcastril, 4 years ago

Add Mixed Precision support by Oriol Tintó

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