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 branches/UKMO/dev_r5518_obs_oper_update_SlaAssim/NEMOGCM/NEMO/OPA_SRC/OBS – NEMO

source: branches/UKMO/dev_r5518_obs_oper_update_SlaAssim/NEMOGCM/NEMO/OPA_SRC/OBS/diaobs.F90 @ 12035

Last change on this file since 12035 was 12035, checked in by kingr, 4 years ago

Added option to remove tides from SLA by taking average over 24h50m.

File size: 95.8 KB
Line 
1MODULE diaobs
2   !!======================================================================
3   !!                       ***  MODULE diaobs  ***
4   !! Observation diagnostics: Computation of the misfit between data and
5   !!                          their model equivalent
6   !!======================================================================
7
8   !!----------------------------------------------------------------------
9   !!   dia_obs_init : Reading and prepare observations
10   !!   dia_obs      : Compute model equivalent to observations
11   !!   dia_obs_wri  : Write observational diagnostics
12   !!   ini_date     : Compute the initial date YYYYMMDD.HHMMSS
13   !!   fin_date     : Compute the final date YYYYMMDD.HHMMSS
14   !!----------------------------------------------------------------------
15   !! * Modules used
16   USE wrk_nemo                 ! Memory Allocation
17   USE par_kind                 ! Precision variables
18   USE in_out_manager           ! I/O manager
19   USE par_oce
20   USE dom_oce                  ! Ocean space and time domain variables
21   USE obs_read_prof            ! Reading and allocation of profile obs
22   USE obs_read_surf            ! Reading and allocation of surface obs
23   USE obs_readmdt              ! Reading and allocation of MDT for SLA.
24   USE obs_prep                 ! Preparation of obs. (grid search etc).
25   USE obs_oper                 ! Observation operators
26   USE obs_write                ! Writing of observation related diagnostics
27   USE obs_grid                 ! Grid searching
28   USE obs_read_altbias         ! Bias treatment for altimeter
29   USE obs_sstbias              ! Bias correction routine for SST
30   USE obs_profiles_def         ! Profile data definitions
31   USE obs_surf_def             ! Surface data definitions
32   USE obs_types                ! Definitions for observation types
33   USE mpp_map                  ! MPP mapping
34   USE lib_mpp                  ! For ctl_warn/stop
35   USE tradmp                   ! For climatological temperature & salinity
36
37   IMPLICIT NONE
38
39   !! * Routine accessibility
40   PRIVATE
41   PUBLIC dia_obs_init, &  ! Initialize and read observations
42      &   dia_obs,      &  ! Compute model equivalent to observations
43      &   dia_obs_wri,  &  ! Write model equivalent to observations
44      &   dia_obs_dealloc  ! Deallocate dia_obs data
45
46   !! * Module variables
47   LOGICAL, PUBLIC :: &
48      &       lk_diaobs = .TRUE.   !: Include this for backwards compatibility at NEMO 3.6.
49   LOGICAL :: ln_diaobs            !: Logical switch for the obs operator
50   LOGICAL :: ln_sstnight          !: Logical switch for night mean SST obs
51   LOGICAL :: ln_default_fp_indegs !: T=> Default obs footprint size specified in degrees, F=> in metres
52   LOGICAL :: ln_sla_fp_indegs     !: T=>     SLA obs footprint size specified in degrees, F=> in metres
53   LOGICAL :: ln_sst_fp_indegs     !: T=>     SST obs footprint size specified in degrees, F=> in metres
54   LOGICAL :: ln_sss_fp_indegs     !: T=>     SSS obs footprint size specified in degrees, F=> in metres
55   LOGICAL :: ln_sic_fp_indegs     !: T=> sea-ice obs footprint size specified in degrees, F=> in metres
56   LOGICAL :: ln_output_clim       !: Logical switch for interpolating and writing T/S climatology
57   LOGICAL :: ln_time_mean_sla_bkg !: Logical switch for applying time mean of SLA background to remove tidal signal
58
59   REAL(wp) :: rn_default_avglamscl !: Default E/W diameter of observation footprint
60   REAL(wp) :: rn_default_avgphiscl !: Default N/S diameter of observation footprint
61   REAL(wp) :: rn_sla_avglamscl     !: E/W diameter of SLA observation footprint
62   REAL(wp) :: rn_sla_avgphiscl     !: N/S diameter of SLA observation footprint
63   REAL(wp) :: rn_sst_avglamscl     !: E/W diameter of SST observation footprint
64   REAL(wp) :: rn_sst_avgphiscl     !: N/S diameter of SST observation footprint
65   REAL(wp) :: rn_sss_avglamscl     !: E/W diameter of SSS observation footprint
66   REAL(wp) :: rn_sss_avgphiscl     !: N/S diameter of SSS observation footprint
67   REAL(wp) :: rn_sic_avglamscl     !: E/W diameter of sea-ice observation footprint
68   REAL(wp) :: rn_sic_avgphiscl     !: N/S diameter of sea-ice observation footprint
69   REAL(wp), PUBLIC :: &
70      &        MeanPeriodHours = 24. + (5./6.) !: Meaning period for surface data.
71
72
73   INTEGER :: nn_1dint         !: Vertical interpolation method
74   INTEGER :: nn_2dint_default !: Default horizontal interpolation method
75   INTEGER :: nn_2dint_sla     !: SLA horizontal interpolation method (-1 = default)
76   INTEGER :: nn_2dint_sst     !: SST horizontal interpolation method (-1 = default)
77   INTEGER :: nn_2dint_sss     !: SSS horizontal interpolation method (-1 = default)
78   INTEGER :: nn_2dint_sic     !: Seaice horizontal interpolation method (-1 = default)
79 
80   INTEGER, DIMENSION(imaxavtypes) :: &
81      & nn_profdavtypes      !: Profile data types representing a daily average
82   INTEGER :: nproftypes     !: Number of profile obs types
83   INTEGER :: nsurftypes     !: Number of surface obs types
84   INTEGER, DIMENSION(:), ALLOCATABLE :: &
85      & nvarsprof, &         !: Number of profile variables
86      & nvarssurf            !: Number of surface variables
87   INTEGER, DIMENSION(:), ALLOCATABLE :: &
88      & nextrprof, &         !: Number of profile extra variables
89      & nextrsurf            !: Number of surface extra variables
90   INTEGER, DIMENSION(:), ALLOCATABLE :: &
91      & n2dintsurf           !: Interpolation option for surface variables
92   REAL(wp), DIMENSION(:), ALLOCATABLE :: &
93      & ravglamscl, &        !: E/W diameter of averaging footprint for surface variables
94      & ravgphiscl           !: N/S diameter of averaging footprint for surface variables
95   LOGICAL, DIMENSION(:), ALLOCATABLE :: &
96      & lfpindegs, &         !: T=> surface obs footprint size specified in degrees, F=> in metres
97      & llnightav            !: Logical for calculating night-time averages
98
99   TYPE(obs_surf), PUBLIC, POINTER, DIMENSION(:) :: &
100      & surfdata, &          !: Initial surface data
101      & surfdataqc           !: Surface data after quality control
102   TYPE(obs_prof), PUBLIC, POINTER, DIMENSION(:) :: &
103      & profdata, &          !: Initial profile data
104      & profdataqc           !: Profile data after quality control
105
106   CHARACTER(len=8), PUBLIC, DIMENSION(:), ALLOCATABLE :: &
107      & cobstypesprof, &     !: Profile obs types
108      & cobstypessurf        !: Surface obs types
109
110   !!----------------------------------------------------------------------
111   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
112   !! $Id$
113   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
114   !!----------------------------------------------------------------------
115
116   !! * Substitutions
117#  include "domzgr_substitute.h90"
118CONTAINS
119
120   SUBROUTINE dia_obs_init
121      !!----------------------------------------------------------------------
122      !!                    ***  ROUTINE dia_obs_init  ***
123      !!         
124      !! ** Purpose : Initialize and read observations
125      !!
126      !! ** Method  : Read the namelist and call reading routines
127      !!
128      !! ** Action  : Read the namelist and call reading routines
129      !!
130      !! History :
131      !!        !  06-03  (K. Mogensen) Original code
132      !!        !  06-05  (A. Weaver) Reformatted
133      !!        !  06-10  (A. Weaver) Cleaning and add controls
134      !!        !  07-03  (K. Mogensen) General handling of profiles
135      !!        !  14-08  (J.While) Incorporated SST bias correction
136      !!        !  15-02  (M. Martin) Simplification of namelist and code
137      !!----------------------------------------------------------------------
138
139      IMPLICIT NONE
140
141      !! * Local declarations
142      INTEGER, PARAMETER :: &
143         & jpmaxnfiles = 1000    ! Maximum number of files for each obs type
144      INTEGER, DIMENSION(:), ALLOCATABLE :: &
145         & ifilesprof, &         ! Number of profile files
146         & ifilessurf            ! Number of surface files
147      INTEGER :: ios             ! Local integer output status for namelist read
148      INTEGER :: jtype           ! Counter for obs types
149      INTEGER :: jvar            ! Counter for variables
150      INTEGER :: jfile           ! Counter for files
151      INTEGER :: jnumsstbias     ! Number of SST bias files to read and apply
152      INTEGER :: n2dint_type     ! Local version of nn_2dint*
153
154      CHARACTER(len=128), DIMENSION(jpmaxnfiles) :: &
155         & cn_profbfiles,      & ! T/S profile input filenames
156         & cn_sstfbfiles,      & ! Sea surface temperature input filenames
157         & cn_slafbfiles,      & ! Sea level anomaly input filenames
158         & cn_sicfbfiles,      & ! Seaice concentration input filenames
159         & cn_velfbfiles,      & ! Velocity profile input filenames
160         & cn_sssfbfiles,      & ! Sea surface salinity input filenames
161         & cn_slchltotfbfiles, & ! Surface total              log10(chlorophyll) input filenames
162         & cn_slchldiafbfiles, & ! Surface diatom             log10(chlorophyll) input filenames
163         & cn_slchlnonfbfiles, & ! Surface non-diatom         log10(chlorophyll) input filenames
164         & cn_slchldinfbfiles, & ! Surface dinoflagellate     log10(chlorophyll) input filenames
165         & cn_slchlmicfbfiles, & ! Surface microphytoplankton log10(chlorophyll) input filenames
166         & cn_slchlnanfbfiles, & ! Surface nanophytoplankton  log10(chlorophyll) input filenames
167         & cn_slchlpicfbfiles, & ! Surface picophytoplankton  log10(chlorophyll) input filenames
168         & cn_schltotfbfiles,  & ! Surface total              chlorophyll        input filenames
169         & cn_slphytotfbfiles, & ! Surface total      log10(phytoplankton carbon) input filenames
170         & cn_slphydiafbfiles, & ! Surface diatom     log10(phytoplankton carbon) input filenames
171         & cn_slphynonfbfiles, & ! Surface non-diatom log10(phytoplankton carbon) input filenames
172         & cn_sspmfbfiles,     & ! Surface suspended particulate matter input filenames
173         & cn_skd490fbfiles,   & ! Surface Kd490 input filenames
174         & cn_sfco2fbfiles,    & ! Surface fugacity         of carbon dioxide input filenames
175         & cn_spco2fbfiles,    & ! Surface partial pressure of carbon dioxide input filenames
176         & cn_plchltotfbfiles, & ! Profile total log10(chlorophyll) input filenames
177         & cn_pchltotfbfiles,  & ! Profile total chlorophyll input filenames
178         & cn_pno3fbfiles,     & ! Profile nitrate input filenames
179         & cn_psi4fbfiles,     & ! Profile silicate input filenames
180         & cn_ppo4fbfiles,     & ! Profile phosphate input filenames
181         & cn_pdicfbfiles,     & ! Profile dissolved inorganic carbon input filenames
182         & cn_palkfbfiles,     & ! Profile alkalinity input filenames
183         & cn_pphfbfiles,      & ! Profile pH input filenames
184         & cn_po2fbfiles,      & ! Profile dissolved oxygen input filenames
185         & cn_sstbiasfiles       ! SST bias input filenames
186
187      CHARACTER(LEN=128) :: &
188         & cn_altbiasfile        ! Altimeter bias input filename
189
190
191      LOGICAL :: ln_t3d          ! Logical switch for temperature profiles
192      LOGICAL :: ln_s3d          ! Logical switch for salinity profiles
193      LOGICAL :: ln_sla          ! Logical switch for sea level anomalies
194      LOGICAL :: ln_sst          ! Logical switch for sea surface temperature
195      LOGICAL :: ln_sic          ! Logical switch for sea ice concentration
196      LOGICAL :: ln_sss          ! Logical switch for sea surface salinity obs
197      LOGICAL :: ln_vel3d        ! Logical switch for velocity (u,v) obs
198      LOGICAL :: ln_slchltot     ! Logical switch for surface total              log10(chlorophyll) obs
199      LOGICAL :: ln_slchldia     ! Logical switch for surface diatom             log10(chlorophyll) obs
200      LOGICAL :: ln_slchlnon     ! Logical switch for surface non-diatom         log10(chlorophyll) obs
201      LOGICAL :: ln_slchldin     ! Logical switch for surface dinoflagellate     log10(chlorophyll) obs
202      LOGICAL :: ln_slchlmic     ! Logical switch for surface microphytoplankton log10(chlorophyll) obs
203      LOGICAL :: ln_slchlnan     ! Logical switch for surface nanophytoplankton  log10(chlorophyll) obs
204      LOGICAL :: ln_slchlpic     ! Logical switch for surface picophytoplankton  log10(chlorophyll) obs
205      LOGICAL :: ln_schltot      ! Logical switch for surface total              chlorophyll        obs
206      LOGICAL :: ln_slphytot     ! Logical switch for surface total      log10(phytoplankton carbon) obs
207      LOGICAL :: ln_slphydia     ! Logical switch for surface diatom     log10(phytoplankton carbon) obs
208      LOGICAL :: ln_slphynon     ! Logical switch for surface non-diatom log10(phytoplankton carbon) obs
209      LOGICAL :: ln_sspm         ! Logical switch for surface suspended particulate matter obs
210      LOGICAL :: ln_skd490       ! Logical switch for surface Kd490
211      LOGICAL :: ln_sfco2        ! Logical switch for surface fugacity         of carbon dioxide obs
212      LOGICAL :: ln_spco2        ! Logical switch for surface partial pressure of carbon dioxide obs
213      LOGICAL :: ln_plchltot     ! Logical switch for profile total log10(chlorophyll) obs
214      LOGICAL :: ln_pchltot      ! Logical switch for profile total chlorophyll obs
215      LOGICAL :: ln_pno3         ! Logical switch for profile nitrate obs
216      LOGICAL :: ln_psi4         ! Logical switch for profile silicate obs
217      LOGICAL :: ln_ppo4         ! Logical switch for profile phosphate obs
218      LOGICAL :: ln_pdic         ! Logical switch for profile dissolved inorganic carbon obs
219      LOGICAL :: ln_palk         ! Logical switch for profile alkalinity obs
220      LOGICAL :: ln_pph          ! Logical switch for profile pH obs
221      LOGICAL :: ln_po2          ! Logical switch for profile dissolved oxygen obs
222      LOGICAL :: ln_nea          ! Logical switch to remove obs near land
223      LOGICAL :: ln_altbias      ! Logical switch for altimeter bias
224      LOGICAL :: ln_sstbias      ! Logical switch for bias correction of SST
225      LOGICAL :: ln_ignmis       ! Logical switch for ignoring missing files
226      LOGICAL :: ln_s_at_t       ! Logical switch to compute model S at T obs
227      LOGICAL :: ln_bound_reject ! Logical switch for rejecting obs near the boundary
228
229      REAL(dp) :: rn_dobsini     ! Obs window start date YYYYMMDD.HHMMSS
230      REAL(dp) :: rn_dobsend     ! Obs window end date   YYYYMMDD.HHMMSS
231
232      REAL(wp) :: ztype_avglamscl ! Local version of rn_*_avglamscl
233      REAL(wp) :: ztype_avgphiscl ! Local version of rn_*_avgphiscl
234
235      CHARACTER(len=128), DIMENSION(:,:), ALLOCATABLE :: &
236         & clproffiles, &        ! Profile filenames
237         & clsurffiles           ! Surface filenames
238      CHARACTER(len=8), DIMENSION(:), ALLOCATABLE :: clvars  ! Expected variable names
239
240      LOGICAL, DIMENSION(:), ALLOCATABLE :: llvar   ! Logical for profile variable read
241      LOGICAL :: ltype_fp_indegs ! Local version of ln_*_fp_indegs
242      LOGICAL :: ltype_night     ! Local version of ln_sstnight (false for other variables)
243      LOGICAL :: ltype_clim      ! Local version of ln_output_clim
244     
245
246      REAL(wp), POINTER, DIMENSION(:,:,:) :: &
247         & zglam                 ! Model longitudes for profile variables
248      REAL(wp), POINTER, DIMENSION(:,:,:) :: &
249         & zgphi                 ! Model latitudes for profile variables
250      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: &
251         & zmask                 ! Model land/sea mask associated with variables
252
253
254      NAMELIST/namobs/ln_diaobs, ln_t3d, ln_s3d, ln_sla,              &
255         &            ln_sst, ln_sic, ln_sss, ln_vel3d,               &
256         &            ln_slchltot, ln_slchldia, ln_slchlnon,          &
257         &            ln_slchldin, ln_slchlmic, ln_slchlnan,          &
258         &            ln_slchlpic, ln_schltot,                        &
259         &            ln_slphytot, ln_slphydia, ln_slphynon,          &
260         &            ln_sspm,     ln_sfco2,    ln_spco2,             &
261         &            ln_skd490,                                      &
262         &            ln_plchltot, ln_pchltot,  ln_pno3,              &
263         &            ln_psi4,     ln_ppo4,     ln_pdic,              &
264         &            ln_palk,     ln_pph,      ln_po2,               &
265         &            ln_altbias, ln_sstbias, ln_nea,                 &
266         &            ln_grid_global, ln_grid_search_lookup,          &
267         &            ln_ignmis, ln_s_at_t, ln_bound_reject,          &
268         &            ln_sstnight,  ln_output_clim,                   &
269         &            ln_time_mean_sla_bkg, ln_default_fp_indegs,     &
270         &            ln_sla_fp_indegs, ln_sst_fp_indegs,             &
271         &            ln_sss_fp_indegs, ln_sic_fp_indegs,             &
272         &            cn_profbfiles, cn_slafbfiles,                   &
273         &            cn_sstfbfiles, cn_sicfbfiles,                   &
274         &            cn_velfbfiles, cn_sssfbfiles,                   &
275         &            cn_slchltotfbfiles, cn_slchldiafbfiles,         &
276         &            cn_slchlnonfbfiles, cn_slchldinfbfiles,         &
277         &            cn_slchlmicfbfiles, cn_slchlnanfbfiles,         &
278         &            cn_slchlpicfbfiles, cn_schltotfbfiles,          &
279         &            cn_slphytotfbfiles, cn_slphydiafbfiles,         &
280         &            cn_slphynonfbfiles, cn_sspmfbfiles,             &
281         &            cn_skd490fbfiles,                               &
282         &            cn_sfco2fbfiles, cn_spco2fbfiles,               &
283         &            cn_plchltotfbfiles, cn_pchltotfbfiles,          &
284         &            cn_pno3fbfiles, cn_psi4fbfiles, cn_ppo4fbfiles, &
285         &            cn_pdicfbfiles, cn_palkfbfiles, cn_pphfbfiles,  &
286         &            cn_po2fbfiles,                                  &
287         &            cn_sstbiasfiles, cn_altbiasfile,                &
288         &            cn_gridsearchfile, rn_gridsearchres,            &
289         &            rn_dobsini, rn_dobsend,                         &
290         &            rn_default_avglamscl, rn_default_avgphiscl,     &
291         &            rn_sla_avglamscl, rn_sla_avgphiscl,             &
292         &            rn_sst_avglamscl, rn_sst_avgphiscl,             &
293         &            rn_sss_avglamscl, rn_sss_avgphiscl,             &
294         &            rn_sic_avglamscl, rn_sic_avgphiscl,             &
295         &            nn_1dint, nn_2dint_default,                     &
296         &            nn_2dint_sla, nn_2dint_sst,                     &
297         &            nn_2dint_sss, nn_2dint_sic,                     &
298         &            nn_msshc, rn_mdtcorr, rn_mdtcutoff,             &
299         &            nn_profdavtypes
300
301      !-----------------------------------------------------------------------
302      ! Read namelist parameters
303      !-----------------------------------------------------------------------
304
305      ! Some namelist arrays need initialising
306      cn_profbfiles(:)      = ''
307      cn_slafbfiles(:)      = ''
308      cn_sstfbfiles(:)      = ''
309      cn_sicfbfiles(:)      = ''
310      cn_velfbfiles(:)      = ''
311      cn_sssfbfiles(:)      = ''
312      cn_slchltotfbfiles(:) = ''
313      cn_slchldiafbfiles(:) = ''
314      cn_slchlnonfbfiles(:) = ''
315      cn_slchldinfbfiles(:) = ''
316      cn_slchlmicfbfiles(:) = ''
317      cn_slchlnanfbfiles(:) = ''
318      cn_slchlpicfbfiles(:) = ''
319      cn_schltotfbfiles(:)  = ''
320      cn_slphytotfbfiles(:) = ''
321      cn_slphydiafbfiles(:) = ''
322      cn_slphynonfbfiles(:) = ''
323      cn_sspmfbfiles(:)     = ''
324      cn_skd490fbfiles(:)   = ''
325      cn_sfco2fbfiles(:)    = ''
326      cn_spco2fbfiles(:)    = ''
327      cn_plchltotfbfiles(:) = ''
328      cn_pchltotfbfiles(:)  = ''
329      cn_pno3fbfiles(:)     = ''
330      cn_psi4fbfiles(:)     = ''
331      cn_ppo4fbfiles(:)     = ''
332      cn_pdicfbfiles(:)     = ''
333      cn_palkfbfiles(:)     = ''
334      cn_pphfbfiles(:)      = ''
335      cn_po2fbfiles(:)      = ''
336      cn_sstbiasfiles(:)    = ''
337      nn_profdavtypes(:)    = -1
338
339      CALL ini_date( rn_dobsini )
340      CALL fin_date( rn_dobsend )
341
342      ! Read namelist namobs : control observation diagnostics
343      REWIND( numnam_ref )   ! Namelist namobs in reference namelist
344      READ  ( numnam_ref, namobs, IOSTAT = ios, ERR = 901)
345901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namobs in reference namelist', lwp )
346
347      REWIND( numnam_cfg )   ! Namelist namobs in configuration namelist
348      READ  ( numnam_cfg, namobs, IOSTAT = ios, ERR = 902 )
349902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namobs in configuration namelist', lwp )
350      IF(lwm) WRITE ( numond, namobs )
351
352      lk_diaobs = .FALSE.
353#if defined key_diaobs
354      IF ( ln_diaobs ) lk_diaobs = .TRUE.
355#endif
356
357      IF ( .NOT. lk_diaobs ) THEN
358         IF(lwp) WRITE(numout,cform_war)
359         IF(lwp) WRITE(numout,*)' ln_diaobs is set to false or key_diaobs is not set, so not calling dia_obs'
360         RETURN
361      ENDIF
362
363      IF(lwp) THEN
364         WRITE(numout,*)
365         WRITE(numout,*) 'dia_obs_init : Observation diagnostic initialization'
366         WRITE(numout,*) '~~~~~~~~~~~~'
367         WRITE(numout,*) '          Namelist namobs : set observation diagnostic parameters' 
368         WRITE(numout,*) '             Logical switch for T profile observations                ln_t3d = ', ln_t3d
369         WRITE(numout,*) '             Logical switch for S profile observations                ln_s3d = ', ln_s3d
370         WRITE(numout,*) '             Logical switch for SLA observations                      ln_sla = ', ln_sla
371         WRITE(numout,*) '             Logical switch for SST observations                      ln_sst = ', ln_sst
372         WRITE(numout,*) '             Logical switch for Sea Ice observations                  ln_sic = ', ln_sic
373         WRITE(numout,*) '             Logical switch for velocity observations               ln_vel3d = ', ln_vel3d
374         WRITE(numout,*) '             Logical switch for SSS observations                      ln_sss = ', ln_sss
375         WRITE(numout,*) '             Logical switch for surface total logchl obs         ln_slchltot = ', ln_slchltot
376         WRITE(numout,*) '             Logical switch for surface diatom logchl obs        ln_slchldia = ', ln_slchldia
377         WRITE(numout,*) '             Logical switch for surface non-diatom logchl obs    ln_slchlnon = ', ln_slchlnon
378         WRITE(numout,*) '             Logical switch for surface dino logchl obs          ln_slchldin = ', ln_slchldin
379         WRITE(numout,*) '             Logical switch for surface micro logchl obs         ln_slchlmic = ', ln_slchlmic
380         WRITE(numout,*) '             Logical switch for surface nano logchl obs          ln_slchlnan = ', ln_slchlnan
381         WRITE(numout,*) '             Logical switch for surface pico logchl obs          ln_slchlpic = ', ln_slchlpic
382         WRITE(numout,*) '             Logical switch for surface total chl obs             ln_schltot = ', ln_schltot
383         WRITE(numout,*) '             Logical switch for surface total log(phyC) obs      ln_slphytot = ', ln_slphytot
384         WRITE(numout,*) '             Logical switch for surface diatom log(phyC) obs     ln_slphydia = ', ln_slphydia
385         WRITE(numout,*) '             Logical switch for surface non-diatom log(phyC) obs ln_slphynon = ', ln_slphynon
386         WRITE(numout,*) '             Logical switch for surface SPM observations             ln_sspm = ', ln_sspm
387         WRITE(numout,*) '             Logical switch for surface Kd490 observations         ln_skd490 = ', ln_skd490
388         WRITE(numout,*) '             Logical switch for surface fCO2 observations           ln_sfco2 = ', ln_sfco2
389         WRITE(numout,*) '             Logical switch for surface pCO2 observations           ln_spco2 = ', ln_spco2
390         WRITE(numout,*) '             Logical switch for profile total logchl obs         ln_plchltot = ', ln_plchltot
391         WRITE(numout,*) '             Logical switch for profile total chl obs             ln_pchltot = ', ln_pchltot
392         WRITE(numout,*) '             Logical switch for profile nitrate obs                  ln_pno3 = ', ln_pno3
393         WRITE(numout,*) '             Logical switch for profile silicate obs                 ln_psi4 = ', ln_psi4
394         WRITE(numout,*) '             Logical switch for profile phosphate obs                ln_ppo4 = ', ln_ppo4
395         WRITE(numout,*) '             Logical switch for profile DIC obs                      ln_pdic = ', ln_pdic
396         WRITE(numout,*) '             Logical switch for profile alkalinity obs               ln_palk = ', ln_palk
397         WRITE(numout,*) '             Logical switch for profile pH obs                        ln_pph = ', ln_pph
398         WRITE(numout,*) '             Logical switch for profile oxygen obs                    ln_po2 = ', ln_po2
399         WRITE(numout,*) '             Global distribution of observations              ln_grid_global = ', ln_grid_global
400         WRITE(numout,*) '             Logical switch for obs grid search lookup ln_grid_search_lookup = ', ln_grid_search_lookup
401         IF (ln_grid_search_lookup) &
402            WRITE(numout,*) '             Grid search lookup file header                cn_gridsearchfile = ', cn_gridsearchfile
403         WRITE(numout,*) '             Initial date in window YYYYMMDD.HHMMSS               rn_dobsini = ', rn_dobsini
404         WRITE(numout,*) '             Final date in window YYYYMMDD.HHMMSS                 rn_dobsend = ', rn_dobsend
405         WRITE(numout,*) '             Type of vertical interpolation method                  nn_1dint = ', nn_1dint
406         WRITE(numout,*) '             Default horizontal interpolation method        nn_2dint_default = ', nn_2dint_default
407         WRITE(numout,*) '             Type of horizontal interpolation method for SLA    nn_2dint_sla = ', nn_2dint_sla
408         WRITE(numout,*) '             Type of horizontal interpolation method for SST    nn_2dint_sst = ', nn_2dint_sst
409         WRITE(numout,*) '             Type of horizontal interpolation method for SSS    nn_2dint_sss = ', nn_2dint_sss
410         WRITE(numout,*) '             Type of horizontal interpolation method for SIC    nn_2dint_sic = ', nn_2dint_sic
411         WRITE(numout,*) '             Default E/W diameter of obs footprint      rn_default_avglamscl = ', rn_default_avglamscl
412         WRITE(numout,*) '             Default N/S diameter of obs footprint      rn_default_avgphiscl = ', rn_default_avgphiscl
413         WRITE(numout,*) '             Default obs footprint in deg [T] or m [F]  ln_default_fp_indegs = ', ln_default_fp_indegs
414         WRITE(numout,*) '             SLA E/W diameter of obs footprint              rn_sla_avglamscl = ', rn_sla_avglamscl
415         WRITE(numout,*) '             SLA N/S diameter of obs footprint              rn_sla_avgphiscl = ', rn_sla_avgphiscl
416         WRITE(numout,*) '             SLA obs footprint in deg [T] or m [F]          ln_sla_fp_indegs = ', ln_sla_fp_indegs
417         WRITE(numout,*) '             SST E/W diameter of obs footprint              rn_sst_avglamscl = ', rn_sst_avglamscl
418         WRITE(numout,*) '             SST N/S diameter of obs footprint              rn_sst_avgphiscl = ', rn_sst_avgphiscl
419         WRITE(numout,*) '             SST obs footprint in deg [T] or m [F]          ln_sst_fp_indegs = ', ln_sst_fp_indegs
420         WRITE(numout,*) '             SIC E/W diameter of obs footprint              rn_sic_avglamscl = ', rn_sic_avglamscl
421         WRITE(numout,*) '             SIC N/S diameter of obs footprint              rn_sic_avgphiscl = ', rn_sic_avgphiscl
422         WRITE(numout,*) '             SIC obs footprint in deg [T] or m [F]          ln_sic_fp_indegs = ', ln_sic_fp_indegs
423         WRITE(numout,*) '             Rejection of observations near land switch               ln_nea = ', ln_nea
424         WRITE(numout,*) '             Rejection of obs near open bdys                 ln_bound_reject = ', ln_bound_reject
425         WRITE(numout,*) '             MSSH correction scheme                                 nn_msshc = ', nn_msshc
426         WRITE(numout,*) '             MDT  correction                                      rn_mdtcorr = ', rn_mdtcorr
427         WRITE(numout,*) '             MDT cutoff for computed correction                 rn_mdtcutoff = ', rn_mdtcutoff
428         WRITE(numout,*) '             Logical switch for alt bias                          ln_altbias = ', ln_altbias
429         WRITE(numout,*) '             Logical switch for sst bias                          ln_sstbias = ', ln_sstbias
430         WRITE(numout,*) '             Logical switch for ignoring missing files             ln_ignmis = ', ln_ignmis
431         WRITE(numout,*) '             Daily average types                             nn_profdavtypes = ', nn_profdavtypes
432         WRITE(numout,*) '             Logical switch for night-time SST obs               ln_sstnight = ', ln_sstnight
433         WRITE(numout,*) '             Logical switch for writing climat. at obs points ln_output_clim = ', ln_output_clim
434         WRITE(numout,*) '             Logical switch for time-mean of SLA        ln_time_mean_sla_bkg = ', ln_time_mean_sla_bkg
435      ENDIF
436      !-----------------------------------------------------------------------
437      ! Set up list of observation types to be used
438      ! and the files associated with each type
439      !-----------------------------------------------------------------------
440
441      nproftypes = COUNT( (/ln_t3d .OR. ln_s3d, ln_vel3d, ln_plchltot,          &
442         &                  ln_pchltot,  ln_pno3,     ln_psi4,     ln_ppo4,     &
443         &                  ln_pdic,     ln_palk,     ln_pph,      ln_po2 /) )
444      nsurftypes = COUNT( (/ln_sla, ln_sst, ln_sic, ln_sss,                     &
445         &                  ln_slchltot, ln_slchldia, ln_slchlnon, ln_slchldin, &
446         &                  ln_slchlmic, ln_slchlnan, ln_slchlpic, ln_schltot,  &
447         &                  ln_slphytot, ln_slphydia, ln_slphynon, ln_sspm,     &
448         &                  ln_skd490,   ln_sfco2,    ln_spco2 /) )
449
450      IF ( nproftypes == 0 .AND. nsurftypes == 0 ) THEN
451         IF(lwp) WRITE(numout,cform_war)
452         IF(lwp) WRITE(numout,*) ' ln_diaobs is set to true, but all obs operator logical flags', &
453            &                    ' are set to .FALSE. so turning off calls to dia_obs'
454         nwarn = nwarn + 1
455         lk_diaobs = .FALSE.
456         RETURN
457      ENDIF
458
459      IF ( ln_output_clim .AND. ( .NOT. ln_tradmp ) ) THEN
460         IF(lwp) WRITE(numout,cform_war)
461         IF(lwp) WRITE(numout,*) ' ln_output_clim is true, but ln_tradmp is false', &
462            &                    ' so climatological T/S not available and will not be output'
463         nwarn = nwarn + 1
464         ln_output_clim = .FALSE.
465      ENDIF
466     
467
468      IF(lwp) WRITE(numout,*) '          Number of profile obs types: ',nproftypes
469      IF ( nproftypes > 0 ) THEN
470
471         ALLOCATE( cobstypesprof(nproftypes) )
472         ALLOCATE( ifilesprof(nproftypes) )
473         ALLOCATE( clproffiles(nproftypes,jpmaxnfiles) )
474
475         jtype = 0
476         IF (ln_t3d .OR. ln_s3d) THEN
477            jtype = jtype + 1
478            cobstypesprof(jtype) = 'prof'
479            clproffiles(jtype,:) = cn_profbfiles
480         ENDIF
481         IF (ln_vel3d) THEN
482            jtype = jtype + 1
483            cobstypesprof(jtype) =  'vel'
484            clproffiles(jtype,:) = cn_velfbfiles
485         ENDIF
486         IF (ln_plchltot) THEN
487            jtype = jtype + 1
488            cobstypesprof(jtype) = 'plchltot'
489            clproffiles(jtype,:) = cn_plchltotfbfiles
490         ENDIF
491         IF (ln_pchltot) THEN
492            jtype = jtype + 1
493            cobstypesprof(jtype) = 'pchltot'
494            clproffiles(jtype,:) = cn_pchltotfbfiles
495         ENDIF
496         IF (ln_pno3) THEN
497            jtype = jtype + 1
498            cobstypesprof(jtype) = 'pno3'
499            clproffiles(jtype,:) = cn_pno3fbfiles
500         ENDIF
501         IF (ln_psi4) THEN
502            jtype = jtype + 1
503            cobstypesprof(jtype) = 'psi4'
504            clproffiles(jtype,:) = cn_psi4fbfiles
505         ENDIF
506         IF (ln_ppo4) THEN
507            jtype = jtype + 1
508            cobstypesprof(jtype) = 'ppo4'
509            clproffiles(jtype,:) = cn_ppo4fbfiles
510         ENDIF
511         IF (ln_pdic) THEN
512            jtype = jtype + 1
513            cobstypesprof(jtype) = 'pdic'
514            clproffiles(jtype,:) = cn_pdicfbfiles
515         ENDIF
516         IF (ln_palk) THEN
517            jtype = jtype + 1
518            cobstypesprof(jtype) = 'palk'
519            clproffiles(jtype,:) = cn_palkfbfiles
520         ENDIF
521         IF (ln_pph) THEN
522            jtype = jtype + 1
523            cobstypesprof(jtype) = 'pph'
524            clproffiles(jtype,:) = cn_pphfbfiles
525         ENDIF
526         IF (ln_po2) THEN
527            jtype = jtype + 1
528            cobstypesprof(jtype) = 'po2'
529            clproffiles(jtype,:) = cn_po2fbfiles
530         ENDIF
531
532         CALL obs_settypefiles( nproftypes, jpmaxnfiles, ifilesprof, cobstypesprof, clproffiles )
533
534      ENDIF
535
536      IF(lwp) WRITE(numout,*)'          Number of surface obs types: ',nsurftypes
537      IF ( nsurftypes > 0 ) THEN
538
539         ALLOCATE( cobstypessurf(nsurftypes) )
540         ALLOCATE( ifilessurf(nsurftypes) )
541         ALLOCATE( clsurffiles(nsurftypes, jpmaxnfiles) )
542         ALLOCATE(n2dintsurf(nsurftypes))
543         ALLOCATE(ravglamscl(nsurftypes))
544         ALLOCATE(ravgphiscl(nsurftypes))
545         ALLOCATE(lfpindegs(nsurftypes))
546         ALLOCATE(llnightav(nsurftypes))
547
548         jtype = 0
549         IF (ln_sla) THEN
550            jtype = jtype + 1
551            cobstypessurf(jtype) = 'sla'
552            clsurffiles(jtype,:) = cn_slafbfiles
553         ENDIF
554         IF (ln_sst) THEN
555            jtype = jtype + 1
556            cobstypessurf(jtype) = 'sst'
557            clsurffiles(jtype,:) = cn_sstfbfiles
558         ENDIF
559         IF (ln_sic) THEN
560            jtype = jtype + 1
561            cobstypessurf(jtype) = 'sic'
562            clsurffiles(jtype,:) = cn_sicfbfiles
563         ENDIF
564         IF (ln_sss) THEN
565            jtype = jtype + 1
566            cobstypessurf(jtype) = 'sss'
567            clsurffiles(jtype,:) = cn_sssfbfiles
568         ENDIF
569         IF (ln_slchltot) THEN
570            jtype = jtype + 1
571            cobstypessurf(jtype) = 'slchltot'
572            clsurffiles(jtype,:) = cn_slchltotfbfiles
573         ENDIF
574         IF (ln_slchldia) THEN
575            jtype = jtype + 1
576            cobstypessurf(jtype) = 'slchldia'
577            clsurffiles(jtype,:) = cn_slchldiafbfiles
578         ENDIF
579         IF (ln_slchlnon) THEN
580            jtype = jtype + 1
581            cobstypessurf(jtype) = 'slchlnon'
582            clsurffiles(jtype,:) = cn_slchlnonfbfiles
583         ENDIF
584         IF (ln_slchldin) THEN
585            jtype = jtype + 1
586            cobstypessurf(jtype) = 'slchldin'
587            clsurffiles(jtype,:) = cn_slchldinfbfiles
588         ENDIF
589         IF (ln_slchlmic) THEN
590            jtype = jtype + 1
591            cobstypessurf(jtype) = 'slchlmic'
592            clsurffiles(jtype,:) = cn_slchlmicfbfiles
593         ENDIF
594         IF (ln_slchlnan) THEN
595            jtype = jtype + 1
596            cobstypessurf(jtype) = 'slchlnan'
597            clsurffiles(jtype,:) = cn_slchlnanfbfiles
598         ENDIF
599         IF (ln_slchlpic) THEN
600            jtype = jtype + 1
601            cobstypessurf(jtype) = 'slchlpic'
602            clsurffiles(jtype,:) = cn_slchlpicfbfiles
603         ENDIF
604         IF (ln_schltot) THEN
605            jtype = jtype + 1
606            cobstypessurf(jtype) = 'schltot'
607            clsurffiles(jtype,:) = cn_schltotfbfiles
608         ENDIF
609         IF (ln_slphytot) THEN
610            jtype = jtype + 1
611            cobstypessurf(jtype) = 'slphytot'
612            clsurffiles(jtype,:) = cn_slphytotfbfiles
613         ENDIF
614         IF (ln_slphydia) THEN
615            jtype = jtype + 1
616            cobstypessurf(jtype) = 'slphydia'
617            clsurffiles(jtype,:) = cn_slphydiafbfiles
618         ENDIF
619         IF (ln_slphynon) THEN
620            jtype = jtype + 1
621            cobstypessurf(jtype) = 'slphynon'
622            clsurffiles(jtype,:) = cn_slphynonfbfiles
623         ENDIF
624         IF (ln_sspm) THEN
625            jtype = jtype + 1
626            cobstypessurf(jtype) = 'sspm'
627            clsurffiles(jtype,:) = cn_sspmfbfiles
628         ENDIF
629         IF (ln_skd490) THEN
630            jtype = jtype + 1
631            cobstypessurf(jtype) = 'skd490'
632            clsurffiles(jtype,:) = cn_skd490fbfiles
633         ENDIF
634         IF (ln_sfco2) THEN
635            jtype = jtype + 1
636            cobstypessurf(jtype) = 'sfco2'
637            clsurffiles(jtype,:) = cn_sfco2fbfiles
638         ENDIF
639         IF (ln_spco2) THEN
640            jtype = jtype + 1
641            cobstypessurf(jtype) = 'spco2'
642            clsurffiles(jtype,:) = cn_spco2fbfiles
643         ENDIF
644
645         CALL obs_settypefiles( nsurftypes, jpmaxnfiles, ifilessurf, cobstypessurf, clsurffiles )
646
647         DO jtype = 1, nsurftypes
648
649            IF ( TRIM(cobstypessurf(jtype)) == 'sla' ) THEN
650               IF ( nn_2dint_sla == -1 ) THEN
651                  n2dint_type  = nn_2dint_default
652               ELSE
653                  n2dint_type  = nn_2dint_sla
654               ENDIF
655               ztype_avglamscl = rn_sla_avglamscl
656               ztype_avgphiscl = rn_sla_avgphiscl
657               ltype_fp_indegs = ln_sla_fp_indegs
658               ltype_night     = .FALSE.
659            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sst' ) THEN
660               IF ( nn_2dint_sst == -1 ) THEN
661                  n2dint_type  = nn_2dint_default
662               ELSE
663                  n2dint_type  = nn_2dint_sst
664               ENDIF
665               ztype_avglamscl = rn_sst_avglamscl
666               ztype_avgphiscl = rn_sst_avgphiscl
667               ltype_fp_indegs = ln_sst_fp_indegs
668               ltype_night     = ln_sstnight
669            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sic' ) THEN
670               IF ( nn_2dint_sic == -1 ) THEN
671                  n2dint_type  = nn_2dint_default
672               ELSE
673                  n2dint_type  = nn_2dint_sic
674               ENDIF
675               ztype_avglamscl = rn_sic_avglamscl
676               ztype_avgphiscl = rn_sic_avgphiscl
677               ltype_fp_indegs = ln_sic_fp_indegs
678               ltype_night     = .FALSE.
679            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sss' ) THEN
680               IF ( nn_2dint_sss == -1 ) THEN
681                  n2dint_type  = nn_2dint_default
682               ELSE
683                  n2dint_type  = nn_2dint_sss
684               ENDIF
685               ztype_avglamscl = rn_sss_avglamscl
686               ztype_avgphiscl = rn_sss_avgphiscl
687               ltype_fp_indegs = ln_sss_fp_indegs
688               ltype_night     = .FALSE.
689            ELSE
690               n2dint_type     = nn_2dint_default
691               ztype_avglamscl = rn_default_avglamscl
692               ztype_avgphiscl = rn_default_avgphiscl
693               ltype_fp_indegs = ln_default_fp_indegs
694               ltype_night     = .FALSE.
695            ENDIF
696           
697            CALL obs_setinterpopts( nsurftypes, jtype, TRIM(cobstypessurf(jtype)), &
698               &                    nn_2dint_default, n2dint_type,                 &
699               &                    ztype_avglamscl, ztype_avgphiscl,              &
700               &                    ltype_fp_indegs, ltype_night,                  &
701               &                    n2dintsurf, ravglamscl, ravgphiscl,            &
702               &                    lfpindegs, llnightav )
703
704         END DO
705
706      ENDIF
707
708      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
709
710
711      !-----------------------------------------------------------------------
712      ! Obs operator parameter checking and initialisations
713      !-----------------------------------------------------------------------
714
715      IF ( ln_vel3d .AND. ( .NOT. ln_grid_global ) ) THEN
716         CALL ctl_stop( 'Velocity data only works with ln_grid_global=.true.' )
717         RETURN
718      ENDIF
719
720      IF ( ( nn_1dint < 0 ) .OR. ( nn_1dint > 1 ) ) THEN
721         CALL ctl_stop(' Choice of vertical (1D) interpolation method', &
722            &                    ' is not available')
723      ENDIF
724
725      IF ( ( nn_2dint_default < 0 ) .OR. ( nn_2dint_default > 6 ) ) THEN
726         CALL ctl_stop(' Choice of default horizontal (2D) interpolation method', &
727            &                    ' is not available')
728      ENDIF
729
730      CALL obs_typ_init
731
732      CALL mppmap_init
733
734      CALL obs_grid_setup( )
735
736      !-----------------------------------------------------------------------
737      ! Depending on switches read the various observation types
738      !-----------------------------------------------------------------------
739
740      IF ( nproftypes > 0 ) THEN
741
742         ALLOCATE(profdata(nproftypes))
743         ALLOCATE(profdataqc(nproftypes))
744         ALLOCATE(nvarsprof(nproftypes))
745         ALLOCATE(nextrprof(nproftypes))
746
747         DO jtype = 1, nproftypes
748           
749            ltype_clim = .FALSE. 
750           
751            IF ( TRIM(cobstypesprof(jtype)) == 'prof' ) THEN
752               nvarsprof(jtype) = 2
753               nextrprof(jtype) = 1
754               IF ( ln_output_clim ) ltype_clim = .TRUE.             
755               ALLOCATE(llvar(nvarsprof(jtype)))
756               ALLOCATE(clvars(nvarsprof(jtype)))
757               CALL wrk_alloc( jpi, jpj,      nvarsprof(jtype), zglam )
758               CALL wrk_alloc( jpi, jpj,      nvarsprof(jtype), zgphi )
759               CALL wrk_alloc( jpi, jpj, jpk, nvarsprof(jtype), zmask )
760               llvar(1)       = ln_t3d
761               llvar(2)       = ln_s3d
762               clvars(1)      = 'POTM'
763               clvars(2)      = 'PSAL'
764               zglam(:,:,1)   = glamt(:,:)
765               zglam(:,:,2)   = glamt(:,:)
766               zgphi(:,:,1)   = gphit(:,:)
767               zgphi(:,:,2)   = gphit(:,:)
768               zmask(:,:,:,1) = tmask(:,:,:)
769               zmask(:,:,:,2) = tmask(:,:,:)
770            ELSE IF ( TRIM(cobstypesprof(jtype)) == 'vel' )  THEN
771               nvarsprof(jtype) = 2
772               nextrprof(jtype) = 2
773               ALLOCATE(llvar(nvarsprof(jtype)))
774               ALLOCATE(clvars(nvarsprof(jtype)))
775               CALL wrk_alloc( jpi, jpj,      nvarsprof(jtype), zglam )
776               CALL wrk_alloc( jpi, jpj,      nvarsprof(jtype), zgphi )
777               CALL wrk_alloc( jpi, jpj, jpk, nvarsprof(jtype), zmask )
778               llvar(1)       = ln_vel3d
779               llvar(2)       = ln_vel3d
780               clvars(1)      = 'UVEL'
781               clvars(2)      = 'VVEL'
782               zglam(:,:,1)   = glamu(:,:)
783               zglam(:,:,2)   = glamv(:,:)
784               zgphi(:,:,1)   = gphiu(:,:)
785               zgphi(:,:,2)   = gphiv(:,:)
786               zmask(:,:,:,1) = umask(:,:,:)
787               zmask(:,:,:,2) = vmask(:,:,:)
788            ELSE
789               nvarsprof(jtype) = 1
790               nextrprof(jtype) = 0
791               ALLOCATE(llvar(nvarsprof(jtype)))
792               ALLOCATE(clvars(nvarsprof(jtype)))
793               CALL wrk_alloc( jpi, jpj,      nvarsprof(jtype), zglam )
794               CALL wrk_alloc( jpi, jpj,      nvarsprof(jtype), zgphi )
795               CALL wrk_alloc( jpi, jpj, jpk, nvarsprof(jtype), zmask )
796               llvar(1)       = .TRUE.
797               zglam(:,:,1)   = glamt(:,:)
798               zgphi(:,:,1)   = gphit(:,:)
799               zmask(:,:,:,1) = tmask(:,:,:)
800               IF ( TRIM(cobstypesprof(jtype)) == 'plchltot' )  THEN
801                  clvars(1) = 'PLCHLTOT'
802               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'pchltot' )  THEN
803                  clvars(1) = 'PCHLTOT'
804               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'pno3' )  THEN
805                  clvars(1) = 'PNO3'
806               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'psi4' )  THEN
807                  clvars(1) = 'PSI4'
808               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'ppo4' )  THEN
809                  clvars(1) = 'PPO4'
810               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'pdic' )  THEN
811                  clvars(1) = 'PDIC'
812               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'palk' )  THEN
813                  clvars(1) = 'PALK'
814               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'pph' )  THEN
815                  clvars(1) = 'PPH'
816               ELSE IF ( TRIM(cobstypesprof(jtype)) == 'po2' )  THEN
817                  clvars(1) = 'PO2'
818               ENDIF
819            ENDIF
820
821            !Read in profile or profile obs types
822            CALL obs_rea_prof( profdata(jtype), ifilesprof(jtype),       &
823               &               clproffiles(jtype,1:ifilesprof(jtype)), &
824               &               nvarsprof(jtype), nextrprof(jtype), nitend-nit000+2, &
825               &               rn_dobsini, rn_dobsend, llvar, &
826               &               ln_ignmis, ln_s_at_t, .FALSE., ltype_clim, clvars, &
827               &               kdailyavtypes = nn_profdavtypes )
828
829            DO jvar = 1, nvarsprof(jtype)
830               CALL obs_prof_staend( profdata(jtype), jvar )
831            END DO
832
833            CALL obs_pre_prof( profdata(jtype), profdataqc(jtype), &
834               &               llvar, &
835               &               jpi, jpj, jpk, &
836               &               zmask, zglam, zgphi,  &
837               &               ln_nea, ln_bound_reject, &
838               &               kdailyavtypes = nn_profdavtypes )
839           
840            DEALLOCATE( llvar, clvars )
841            CALL wrk_dealloc( jpi, jpj,      nvarsprof(jtype), zglam )
842            CALL wrk_dealloc( jpi, jpj,      nvarsprof(jtype), zgphi )
843            CALL wrk_dealloc( jpi, jpj, jpk, nvarsprof(jtype), zmask )
844
845         END DO
846
847         DEALLOCATE( ifilesprof, clproffiles )
848
849      ENDIF
850
851      IF ( nsurftypes > 0 ) THEN
852
853         ALLOCATE(surfdata(nsurftypes))
854         ALLOCATE(surfdataqc(nsurftypes))
855         ALLOCATE(nvarssurf(nsurftypes))
856         ALLOCATE(nextrsurf(nsurftypes))
857
858         DO jtype = 1, nsurftypes
859
860            ltype_clim = .FALSE.
861            IF ( ln_output_clim .AND. &
862               & ( ( TRIM(cobstypessurf(jtype)) == 'sst' ) .OR. &
863               &   ( TRIM(cobstypessurf(jtype)) == 'sss' ) ) ) &
864               & ltype_clim = .TRUE.
865
866            IF ( TRIM(cobstypessurf(jtype)) == 'sla' ) THEN
867               nvarssurf(jtype) = 1
868               nextrsurf(jtype) = 2
869            ELSE
870               nvarssurf(jtype) = 1
871               nextrsurf(jtype) = 0
872            ENDIF
873           
874            ALLOCATE( clvars( nvarssurf(jtype) ) )
875
876            IF ( TRIM(cobstypessurf(jtype)) == 'sla' ) THEN
877               clvars(1) = 'SLA'
878            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sst' ) THEN
879               clvars(1) = 'SST'
880            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sic' ) THEN
881               clvars(1) = 'ICECONC'
882            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sss' ) THEN
883               clvars(1) = 'SSS'
884            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slchltot' ) THEN
885               clvars(1) = 'SLCHLTOT'
886            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slchldia' ) THEN
887               clvars(1) = 'SLCHLDIA'
888            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slchlnon' ) THEN
889               clvars(1) = 'SLCHLNON'
890            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slchldin' ) THEN
891               clvars(1) = 'SLCHLDIN'
892            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slchlmic' ) THEN
893               clvars(1) = 'SLCHLMIC'
894            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slchlnan' ) THEN
895               clvars(1) = 'SLCHLNAN'
896            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slchlpic' ) THEN
897               clvars(1) = 'SLCHLPIC'
898            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'schltot' ) THEN
899               clvars(1) = 'SCHLTOT'
900            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slphytot' ) THEN
901               clvars(1) = 'SLPHYTOT'
902            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slphydia' ) THEN
903               clvars(1) = 'SLPHYDIA'
904            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'slphynon' ) THEN
905               clvars(1) = 'SLPHYNON'
906            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sspm' ) THEN
907               clvars(1) = 'SSPM'
908            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'skd490' ) THEN
909               clvars(1) = 'SKD490'
910            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'sfco2' ) THEN
911               clvars(1) = 'SFCO2'
912            ELSE IF ( TRIM(cobstypessurf(jtype)) == 'spco2' ) THEN
913               clvars(1) = 'SPCO2'
914            ENDIF
915
916            !Read in surface obs types
917            CALL obs_rea_surf( surfdata(jtype), ifilessurf(jtype), &
918               &               clsurffiles(jtype,1:ifilessurf(jtype)), &
919               &               nvarssurf(jtype), nextrsurf(jtype), nitend-nit000+2, &
920               &               rn_dobsini, rn_dobsend, MeanPeriodHours, ln_ignmis, .FALSE., &
921               &               llnightav(jtype), ltype_clim, clvars )
922
923            CALL obs_pre_surf( surfdata(jtype), surfdataqc(jtype), ln_nea, ln_bound_reject )
924
925            IF ( TRIM(cobstypessurf(jtype)) == 'sla' ) THEN
926               CALL obs_rea_mdt( surfdataqc(jtype), n2dintsurf(jtype) )
927               IF ( ln_altbias ) &
928                  & CALL obs_rea_altbias ( surfdataqc(jtype), n2dintsurf(jtype), cn_altbiasfile )
929            ENDIF
930
931            IF ( TRIM(cobstypessurf(jtype)) == 'sst' .AND. ln_sstbias ) THEN
932               jnumsstbias = 0
933               DO jfile = 1, jpmaxnfiles
934                  IF ( TRIM(cn_sstbiasfiles(jfile)) /= '' ) &
935                     &  jnumsstbias = jnumsstbias + 1
936               END DO
937               IF ( jnumsstbias == 0 ) THEN
938                  CALL ctl_stop("ln_sstbias set but no bias files to read in")   
939               ENDIF
940
941               CALL obs_app_sstbias( surfdataqc(jtype), n2dintsurf(jtype), & 
942                  &                  jnumsstbias, cn_sstbiasfiles(1:jnumsstbias) ) 
943
944            ENDIF
945           
946            DEALLOCATE( clvars )
947
948         END DO
949
950         DEALLOCATE( ifilessurf, clsurffiles )
951
952      ENDIF
953
954   END SUBROUTINE dia_obs_init
955
956   SUBROUTINE dia_obs( kstp )
957      !!----------------------------------------------------------------------
958      !!                    ***  ROUTINE dia_obs  ***
959      !!         
960      !! ** Purpose : Call the observation operators on each time step
961      !!
962      !! ** Method  : Call the observation operators on each time step to
963      !!              compute the model equivalent of the following data:
964      !!               - Profile data, currently T/S or U/V
965      !!               - Surface data, currently SST, SLA or sea-ice concentration.
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      !!        !  07-03  (K. Mogensen) General handling of profiles
974      !!        !  07-04  (G. Smith) Generalized surface operators
975      !!        !  08-10  (M. Valdivieso) obs operator for velocity profiles
976      !!        !  15-08  (M. Martin) Combined surface/profile routines.
977      !!----------------------------------------------------------------------
978      !! * Modules used
979      USE phycst, ONLY : &         ! Physical constants
980#if defined key_fabm
981         & rt0,          &
982#endif
983         & rday
984      USE oce, ONLY : &            ! Ocean dynamics and tracers variables
985         & tsn,       &
986         & un,        &
987         & vn,        &
988         & sshn
989#if defined  key_lim3
990      USE ice, ONLY : &            ! LIM3 Ice model variables
991         & frld
992#endif
993#if defined key_lim2
994      USE ice_2, ONLY : &          ! LIM2 Ice model variables
995         & frld
996#endif
997#if defined key_cice
998      USE sbc_oce, ONLY : fr_i     ! ice fraction
999#endif
1000#if defined key_top
1001      USE trc, ONLY :  &           ! Biogeochemical state variables
1002         & trn
1003#endif
1004#if defined key_hadocc
1005      USE par_hadocc               ! HadOCC parameters
1006      USE trc, ONLY :  &
1007         & HADOCC_CHL, &
1008         & HADOCC_FCO2, &
1009         & HADOCC_PCO2, &
1010         & HADOCC_FILL_FLT
1011      USE had_bgc_const, ONLY: c2n_p
1012#elif defined key_medusa
1013      USE par_medusa               ! MEDUSA parameters
1014      USE sms_medusa, ONLY: &
1015         & xthetapn, &
1016         & xthetapd
1017#if defined key_roam
1018      USE sms_medusa, ONLY: &
1019         & f2_pco2w, &
1020         & f2_fco2w, &
1021         & f3_pH
1022#endif
1023#elif defined key_fabm
1024      USE par_fabm                 ! FABM parameters
1025      USE fabm, ONLY: &
1026         & fabm_get_interior_diagnostic_data
1027#endif
1028#if defined key_spm
1029      USE par_spm, ONLY: &         ! Sediment parameters
1030         & jp_spm
1031#endif
1032
1033      IMPLICIT NONE
1034
1035      !! * Arguments
1036      INTEGER, INTENT(IN) :: kstp  ! Current timestep
1037      !! * Local declarations
1038      INTEGER :: idaystp           ! Number of timesteps per day
1039      INTEGER :: imeanstp          ! Number of timesteps for sla averaging
1040      INTEGER :: jtype             ! Data loop variable
1041      INTEGER :: jvar              ! Variable number
1042      INTEGER :: ji, jj, jk        ! Loop counters
1043      REAL(wp) :: tiny             ! small number
1044      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: &
1045         & zprofvar, &             ! Model values for variables in a prof ob
1046         & zprofclim               ! Climatology values for variables in a prof ob
1047      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: &
1048         & zprofmask               ! Mask associated with zprofvar
1049      REAL(wp), POINTER, DIMENSION(:,:) :: &
1050         & zsurfvar, &             ! Model values equivalent to surface ob.
1051         & zsurfclim, &            ! Climatology values for variables in a surface ob.
1052         & zsurfmask               ! Mask associated with surface variable
1053      REAL(wp), POINTER, DIMENSION(:,:,:) :: &
1054         & zglam,    &             ! Model longitudes for prof variables
1055         & zgphi                   ! Model latitudes for prof variables
1056      LOGICAL :: llog10            ! Perform log10 transform of variable
1057#if defined key_fabm
1058      REAL(wp), POINTER, DIMENSION(:,:,:) :: &
1059         & fabm_3d                 ! 3D variable from FABM
1060#endif
1061     
1062      IF(lwp) THEN
1063         WRITE(numout,*)
1064         WRITE(numout,*) 'dia_obs : Call the observation operators', kstp
1065         WRITE(numout,*) '~~~~~~~'
1066         CALL FLUSH(numout)
1067      ENDIF
1068
1069      idaystp = NINT( rday / rdt )
1070
1071      !-----------------------------------------------------------------------
1072      ! Call the profile and surface observation operators
1073      !-----------------------------------------------------------------------
1074
1075      IF ( nproftypes > 0 ) THEN
1076
1077         DO jtype = 1, nproftypes
1078
1079            ! Allocate local work arrays
1080            CALL wrk_alloc( jpi, jpj, jpk, profdataqc(jtype)%nvar, zprofvar  )
1081            CALL wrk_alloc( jpi, jpj, jpk, profdataqc(jtype)%nvar, zprofmask )
1082            CALL wrk_alloc( jpi, jpj,      profdataqc(jtype)%nvar, zglam     )
1083            CALL wrk_alloc( jpi, jpj,      profdataqc(jtype)%nvar, zgphi     )
1084            CALL wrk_alloc( jpi, jpj, jpk, profdataqc(jtype)%nvar, zprofclim )   
1085                             
1086            ! Defaults which might change
1087            DO jvar = 1, profdataqc(jtype)%nvar
1088               zprofmask(:,:,:,jvar) = tmask(:,:,:)
1089               zglam(:,:,jvar)       = glamt(:,:)
1090               zgphi(:,:,jvar)       = gphit(:,:)
1091               zprofclim(:,:,:,jvar) = 0._wp
1092            END DO
1093
1094            SELECT CASE ( TRIM(cobstypesprof(jtype)) )
1095
1096            CASE('prof')
1097               zprofvar(:,:,:,1) = tsn(:,:,:,jp_tem)
1098               zprofvar(:,:,:,2) = tsn(:,:,:,jp_sal)
1099               IF ( ln_output_clim ) THEN         
1100                  zprofclim(:,:,:,1) = tclim(:,:,:)
1101                  zprofclim(:,:,:,2) = sclim(:,:,:)
1102               ENDIF
1103               
1104            CASE('vel')
1105               zprofvar(:,:,:,1) = un(:,:,:)
1106               zprofvar(:,:,:,2) = vn(:,:,:)
1107               zprofmask(:,:,:,1) = umask(:,:,:)
1108               zprofmask(:,:,:,2) = vmask(:,:,:)
1109               zglam(:,:,1) = glamu(:,:)
1110               zglam(:,:,2) = glamv(:,:)
1111               zgphi(:,:,1) = gphiu(:,:)
1112               zgphi(:,:,2) = gphiv(:,:)
1113
1114            CASE('plchltot')
1115#if defined key_hadocc
1116               ! Chlorophyll from HadOCC
1117               zprofvar(:,:,:,1) = HADOCC_CHL(:,:,:)
1118#elif defined key_medusa
1119               ! Add non-diatom and diatom chlorophyll from MEDUSA
1120               zprofvar(:,:,:,1) = trn(:,:,:,jpchn) + trn(:,:,:,jpchd)
1121#elif defined key_fabm
1122               ! Add all chlorophyll groups from ERSEM
1123               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_chl1) + trn(:,:,:,jp_fabm_m1+jp_fabm_chl2) + &
1124                  &                trn(:,:,:,jp_fabm_m1+jp_fabm_chl3) + trn(:,:,:,jp_fabm_m1+jp_fabm_chl4)
1125#else
1126               CALL ctl_stop( ' Trying to run plchltot observation operator', &
1127                  &           ' but no biogeochemical model appears to have been defined' )
1128#endif
1129               ! Take the log10 where we can, otherwise exclude
1130               tiny = 1.0e-20
1131               WHERE(zprofvar(:,:,:,:) > tiny .AND. zprofvar(:,:,:,:) /= obfillflt )
1132                  zprofvar(:,:,:,:)  = LOG10(zprofvar(:,:,:,:))
1133               ELSEWHERE
1134                  zprofvar(:,:,:,:)  = obfillflt
1135                  zprofmask(:,:,:,:) = 0
1136               END WHERE
1137               ! Mask out model below any excluded values,
1138               ! to avoid interpolation issues
1139               DO jvar = 1, profdataqc(jtype)%nvar
1140                 DO jj = 1, jpj
1141                    DO ji = 1, jpi
1142                       depth_loop: DO jk = 1, jpk
1143                          IF ( zprofmask(ji,jj,jk,jvar) == 0 ) THEN
1144                             zprofmask(ji,jj,jk:jpk,jvar) = 0
1145                             EXIT depth_loop
1146                          ENDIF
1147                       END DO depth_loop
1148                    END DO
1149                 END DO
1150              END DO
1151
1152            CASE('pchltot')
1153#if defined key_hadocc
1154               ! Chlorophyll from HadOCC
1155               zprofvar(:,:,:,1) = HADOCC_CHL(:,:,:)
1156#elif defined key_medusa
1157               ! Add non-diatom and diatom chlorophyll from MEDUSA
1158               zprofvar(:,:,:,1) = trn(:,:,:,jpchn) + trn(:,:,:,jpchd)
1159#elif defined key_fabm
1160               ! Add all chlorophyll groups from ERSEM
1161               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_chl1) + trn(:,:,:,jp_fabm_m1+jp_fabm_chl2) + &
1162                  &                trn(:,:,:,jp_fabm_m1+jp_fabm_chl3) + trn(:,:,:,jp_fabm_m1+jp_fabm_chl4)
1163#else
1164               CALL ctl_stop( ' Trying to run pchltot observation operator', &
1165                  &           ' but no biogeochemical model appears to have been defined' )
1166#endif
1167
1168            CASE('pno3')
1169#if defined key_hadocc
1170               ! Dissolved inorganic nitrogen from HadOCC
1171               zprofvar(:,:,:,1) = trn(:,:,:,jp_had_nut)
1172#elif defined key_medusa
1173               ! Dissolved inorganic nitrogen from MEDUSA
1174               zprofvar(:,:,:,1) = trn(:,:,:,jpdin)
1175#elif defined key_fabm
1176               ! Nitrate from ERSEM
1177               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_n3n)
1178#else
1179               CALL ctl_stop( ' Trying to run pno3 observation operator', &
1180                  &           ' but no biogeochemical model appears to have been defined' )
1181#endif
1182
1183            CASE('psi4')
1184#if defined key_hadocc
1185               CALL ctl_stop( ' Trying to run psi4 observation operator', &
1186                  &           ' but HadOCC does not simulate silicate' )
1187#elif defined key_medusa
1188               ! Silicate from MEDUSA
1189               zprofvar(:,:,:,1) = trn(:,:,:,jpsil)
1190#elif defined key_fabm
1191               ! Silicate from ERSEM
1192               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_n5s)
1193#else
1194               CALL ctl_stop( ' Trying to run psi4 observation operator', &
1195                  &           ' but no biogeochemical model appears to have been defined' )
1196#endif
1197
1198            CASE('ppo4')
1199#if defined key_hadocc
1200               CALL ctl_stop( ' Trying to run ppo4 observation operator', &
1201                  &           ' but HadOCC does not simulate phosphate' )
1202#elif defined key_medusa
1203               CALL ctl_stop( ' Trying to run ppo4 observation operator', &
1204                  &           ' but MEDUSA does not simulate phosphate' )
1205#elif defined key_fabm
1206               ! Phosphate from ERSEM
1207               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_n1p)
1208#else
1209               CALL ctl_stop( ' Trying to run ppo4 observation operator', &
1210                  &           ' but no biogeochemical model appears to have been defined' )
1211#endif
1212
1213            CASE('pdic')
1214#if defined key_hadocc
1215               ! Dissolved inorganic carbon from HadOCC
1216               zprofvar(:,:,:,1) = trn(:,:,:,jp_had_dic)
1217#elif defined key_medusa
1218               ! Dissolved inorganic carbon from MEDUSA
1219               zprofvar(:,:,:,1) = trn(:,:,:,jpdic)
1220#elif defined key_fabm
1221               ! Dissolved inorganic carbon from ERSEM
1222               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_o3c)
1223#else
1224               CALL ctl_stop( ' Trying to run pdic observation operator', &
1225                  &           ' but no biogeochemical model appears to have been defined' )
1226#endif
1227
1228            CASE('palk')
1229#if defined key_hadocc
1230               ! Alkalinity from HadOCC
1231               zprofvar(:,:,:,1) = trn(:,:,:,jp_had_alk)
1232#elif defined key_medusa
1233               ! Alkalinity from MEDUSA
1234               zprofvar(:,:,:,1) = trn(:,:,:,jpalk)
1235#elif defined key_fabm
1236               ! Alkalinity from ERSEM
1237               zprofvar(:,:,:,1) = fabm_get_interior_diagnostic_data(model, jp_fabm_o3ta)
1238#else
1239               CALL ctl_stop( ' Trying to run palk observation operator', &
1240                  &           ' but no biogeochemical model appears to have been defined' )
1241#endif
1242
1243            CASE('pph')
1244#if defined key_hadocc
1245               CALL ctl_stop( ' Trying to run pph observation operator', &
1246                  &           ' but HadOCC has no pH diagnostic defined' )
1247#elif defined key_medusa && defined key_roam
1248               ! pH from MEDUSA
1249               zprofvar(:,:,:,1) = f3_pH(:,:,:)
1250#elif defined key_fabm
1251               ! pH from ERSEM
1252               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_o3ph)
1253#else
1254               CALL ctl_stop( ' Trying to run pph observation operator', &
1255                  &           ' but no biogeochemical model appears to have been defined' )
1256#endif
1257
1258            CASE('po2')
1259#if defined key_hadocc
1260               CALL ctl_stop( ' Trying to run po2 observation operator', &
1261                  &           ' but HadOCC does not simulate oxygen' )
1262#elif defined key_medusa
1263               ! Oxygen from MEDUSA
1264               zprofvar(:,:,:,1) = trn(:,:,:,jpoxy)
1265#elif defined key_fabm
1266               ! Oxygen from ERSEM
1267               zprofvar(:,:,:,1) = trn(:,:,:,jp_fabm_m1+jp_fabm_o2o)
1268#else
1269               CALL ctl_stop( ' Trying to run po2 observation operator', &
1270                  &           ' but no biogeochemical model appears to have been defined' )
1271#endif
1272
1273            CASE DEFAULT
1274               CALL ctl_stop( 'Unknown profile observation type '//TRIM(cobstypesprof(jtype))//' in dia_obs' )
1275
1276            END SELECT
1277
1278            DO jvar = 1, profdataqc(jtype)%nvar
1279               CALL obs_prof_opt( profdataqc(jtype), kstp, jpi, jpj, jpk,  &
1280                  &               nit000, idaystp, jvar,                   &
1281                  &               zprofvar(:,:,:,jvar),                    &
1282                  &               zprofclim(:,:,:,jvar),                   &
1283                  &               fsdept(:,:,:), fsdepw(:,:,:),            & 
1284                  &               zprofmask(:,:,:,jvar),                   &
1285                  &               zglam(:,:,jvar), zgphi(:,:,jvar),        &
1286                  &               nn_1dint, nn_2dint_default,              &
1287                  &               kdailyavtypes = nn_profdavtypes )
1288            END DO
1289
1290            CALL wrk_dealloc( jpi, jpj, jpk, profdataqc(jtype)%nvar, zprofvar  )
1291            CALL wrk_dealloc( jpi, jpj, jpk, profdataqc(jtype)%nvar, zprofmask )
1292            CALL wrk_dealloc( jpi, jpj,      profdataqc(jtype)%nvar, zglam     )
1293            CALL wrk_dealloc( jpi, jpj,      profdataqc(jtype)%nvar, zgphi     )
1294            CALL wrk_dealloc( jpi, jpj, jpk, profdataqc(jtype)%nvar, zprofclim  )           
1295
1296         END DO
1297
1298      ENDIF
1299
1300      IF ( nsurftypes > 0 ) THEN
1301
1302         !Allocate local work arrays
1303         CALL wrk_alloc( jpi, jpj, zsurfvar )
1304         CALL wrk_alloc( jpi, jpj, zsurfclim )         
1305         CALL wrk_alloc( jpi, jpj, zsurfmask )
1306#if defined key_fabm
1307         CALL wrk_alloc( jpi, jpj, jpk, fabm_3d )
1308#endif
1309
1310         DO jtype = 1, nsurftypes
1311
1312            !Defaults which might be changed
1313            zsurfmask(:,:) = tmask(:,:,1)
1314            zsurfclim(:,:) = 0._wp         
1315            llog10 = .FALSE.
1316
1317            SELECT CASE ( TRIM(cobstypessurf(jtype)) )
1318            CASE('sst')
1319               zsurfvar(:,:) = tsn(:,:,1,jp_tem)
1320               IF ( ln_output_clim ) zsurfclim(:,:) = tclim(:,:,1)
1321            CASE('sla')
1322               zsurfvar(:,:) = sshn(:,:)
1323            CASE('sss')
1324               zsurfvar(:,:) = tsn(:,:,1,jp_sal)
1325               IF ( ln_output_clim ) zsurfclim(:,:) = sclim(:,:,1)             
1326            CASE('sic')
1327               IF ( kstp == 0 ) THEN
1328                  IF ( lwp .AND. surfdataqc(jtype)%nsstpmpp(1) > 0 ) THEN
1329                     CALL ctl_warn( 'Sea-ice not initialised on zeroth '// &
1330                        &           'time-step but some obs are valid then.' )
1331                     WRITE(numout,*)surfdataqc(jtype)%nsstpmpp(1), &
1332                        &           ' sea-ice obs will be missed'
1333                  ENDIF
1334                  surfdataqc(jtype)%nsurfup = surfdataqc(jtype)%nsurfup + &
1335                     &                        surfdataqc(jtype)%nsstp(1)
1336                  CYCLE
1337               ELSE
1338#if defined key_cice
1339                  zsurfvar(:,:) = fr_i(:,:)
1340#elif defined key_lim2 || defined key_lim3
1341                  zsurfvar(:,:) = 1._wp - frld(:,:)
1342#else
1343               CALL ctl_stop( ' Trying to run sea-ice observation operator', &
1344                  &           ' but no sea-ice model appears to have been defined' )
1345#endif
1346               ENDIF
1347
1348            CASE('slchltot')
1349#if defined key_hadocc
1350               ! Surface chlorophyll from HadOCC
1351               zsurfvar(:,:) = HADOCC_CHL(:,:,1)
1352#elif defined key_medusa
1353               ! Add non-diatom and diatom surface chlorophyll from MEDUSA
1354               zsurfvar(:,:) = trn(:,:,1,jpchn) + trn(:,:,1,jpchd)
1355#elif defined key_fabm
1356               ! Add all surface chlorophyll groups from ERSEM
1357               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl1) + trn(:,:,1,jp_fabm_m1+jp_fabm_chl2) + &
1358                  &            trn(:,:,1,jp_fabm_m1+jp_fabm_chl3) + trn(:,:,1,jp_fabm_m1+jp_fabm_chl4)
1359#else
1360               CALL ctl_stop( ' Trying to run slchltot observation operator', &
1361                  &           ' but no biogeochemical model appears to have been defined' )
1362#endif
1363               llog10 = .TRUE.
1364
1365            CASE('slchldia')
1366#if defined key_hadocc
1367               CALL ctl_stop( ' Trying to run slchldia observation operator', &
1368                  &           ' but HadOCC does not explicitly simulate diatoms' )
1369#elif defined key_medusa
1370               ! Diatom surface chlorophyll from MEDUSA
1371               zsurfvar(:,:) = trn(:,:,1,jpchd)
1372#elif defined key_fabm
1373               ! Diatom surface chlorophyll from ERSEM
1374               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl1)
1375#else
1376               CALL ctl_stop( ' Trying to run slchldia observation operator', &
1377                  &           ' but no biogeochemical model appears to have been defined' )
1378#endif
1379               llog10 = .TRUE.
1380
1381            CASE('slchlnon')
1382#if defined key_hadocc
1383               CALL ctl_stop( ' Trying to run slchlnon observation operator', &
1384                  &           ' but HadOCC does not explicitly simulate non-diatoms' )
1385#elif defined key_medusa
1386               ! Non-diatom surface chlorophyll from MEDUSA
1387               zsurfvar(:,:) = trn(:,:,1,jpchn)
1388#elif defined key_fabm
1389               ! Add all non-diatom surface chlorophyll groups from ERSEM
1390               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl2) + &
1391                  &            trn(:,:,1,jp_fabm_m1+jp_fabm_chl3) + trn(:,:,1,jp_fabm_m1+jp_fabm_chl4)
1392#else
1393               CALL ctl_stop( ' Trying to run slchlnon observation operator', &
1394                  &           ' but no biogeochemical model appears to have been defined' )
1395#endif
1396               llog10 = .TRUE.
1397
1398            CASE('slchldin')
1399#if defined key_hadocc
1400               CALL ctl_stop( ' Trying to run slchldin observation operator', &
1401                  &           ' but HadOCC does not explicitly simulate dinoflagellates' )
1402#elif defined key_medusa
1403               CALL ctl_stop( ' Trying to run slchldin observation operator', &
1404                  &           ' but MEDUSA does not explicitly simulate dinoflagellates' )
1405#elif defined key_fabm
1406               ! Dinoflagellate surface chlorophyll from ERSEM
1407               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl4)
1408#else
1409               CALL ctl_stop( ' Trying to run slchldin observation operator', &
1410                  &           ' but no biogeochemical model appears to have been defined' )
1411#endif
1412               llog10 = .TRUE.
1413
1414            CASE('slchlmic')
1415#if defined key_hadocc
1416               CALL ctl_stop( ' Trying to run slchlmic observation operator', &
1417                  &           ' but HadOCC does not explicitly simulate microphytoplankton' )
1418#elif defined key_medusa
1419               CALL ctl_stop( ' Trying to run slchlmic observation operator', &
1420                  &           ' but MEDUSA does not explicitly simulate microphytoplankton' )
1421#elif defined key_fabm
1422               ! Add diatom and dinoflagellate surface chlorophyll from ERSEM
1423               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl1) + trn(:,:,1,jp_fabm_m1+jp_fabm_chl4)
1424#else
1425               CALL ctl_stop( ' Trying to run slchlmic observation operator', &
1426                  &           ' but no biogeochemical model appears to have been defined' )
1427#endif
1428               llog10 = .TRUE.
1429
1430            CASE('slchlnan')
1431#if defined key_hadocc
1432               CALL ctl_stop( ' Trying to run slchlnan observation operator', &
1433                  &           ' but HadOCC does not explicitly simulate nanophytoplankton' )
1434#elif defined key_medusa
1435               CALL ctl_stop( ' Trying to run slchlnan observation operator', &
1436                  &           ' but MEDUSA does not explicitly simulate nanophytoplankton' )
1437#elif defined key_fabm
1438               ! Nanophytoplankton surface chlorophyll from ERSEM
1439               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl2)
1440#else
1441               CALL ctl_stop( ' Trying to run slchlnan observation operator', &
1442                  &           ' but no biogeochemical model appears to have been defined' )
1443#endif
1444               llog10 = .TRUE.
1445
1446            CASE('slchlpic')
1447#if defined key_hadocc
1448               CALL ctl_stop( ' Trying to run slchlpic observation operator', &
1449                  &           ' but HadOCC does not explicitly simulate picophytoplankton' )
1450#elif defined key_medusa
1451               CALL ctl_stop( ' Trying to run slchlpic observation operator', &
1452                  &           ' but MEDUSA does not explicitly simulate picophytoplankton' )
1453#elif defined key_fabm
1454               ! Picophytoplankton surface chlorophyll from ERSEM
1455               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl3)
1456#else
1457               CALL ctl_stop( ' Trying to run slchlpic observation operator', &
1458                  &           ' but no biogeochemical model appears to have been defined' )
1459#endif
1460               llog10 = .TRUE.
1461
1462            CASE('schltot')
1463#if defined key_hadocc
1464               ! Surface chlorophyll from HadOCC
1465               zsurfvar(:,:) = HADOCC_CHL(:,:,1)
1466#elif defined key_medusa
1467               ! Add non-diatom and diatom surface chlorophyll from MEDUSA
1468               zsurfvar(:,:) = trn(:,:,1,jpchn) + trn(:,:,1,jpchd)
1469#elif defined key_fabm
1470               ! Add all surface chlorophyll groups from ERSEM
1471               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_chl1) + trn(:,:,1,jp_fabm_m1+jp_fabm_chl2) + &
1472                  &            trn(:,:,1,jp_fabm_m1+jp_fabm_chl3) + trn(:,:,1,jp_fabm_m1+jp_fabm_chl4)
1473#else
1474               CALL ctl_stop( ' Trying to run schltot observation operator', &
1475                  &           ' but no biogeochemical model appears to have been defined' )
1476#endif
1477
1478            CASE('slphytot')
1479#if defined key_hadocc
1480               ! Surface phytoplankton nitrogen from HadOCC multiplied by C:N ratio
1481               zsurfvar(:,:) = trn(:,:,1,jp_had_phy) * c2n_p
1482#elif defined key_medusa
1483               ! Add non-diatom and diatom surface phytoplankton nitrogen from MEDUSA
1484               ! multiplied by C:N ratio for each
1485               zsurfvar(:,:) = (trn(:,:,1,jpphn) * xthetapn) + (trn(:,:,1,jpphd) * xthetapd)
1486#elif defined key_fabm
1487               ! Add all surface phytoplankton carbon groups from ERSEM
1488               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_p1c) + trn(:,:,1,jp_fabm_m1+jp_fabm_p2c) + &
1489                  &            trn(:,:,1,jp_fabm_m1+jp_fabm_p3c) + trn(:,:,1,jp_fabm_m1+jp_fabm_p4c)
1490#else
1491               CALL ctl_stop( ' Trying to run slphytot observation operator', &
1492                  &           ' but no biogeochemical model appears to have been defined' )
1493#endif
1494               llog10 = .TRUE.
1495
1496            CASE('slphydia')
1497#if defined key_hadocc
1498               CALL ctl_stop( ' Trying to run slphydia observation operator', &
1499                  &           ' but HadOCC does not explicitly simulate diatoms' )
1500#elif defined key_medusa
1501               ! Diatom surface phytoplankton nitrogen from MEDUSA multiplied by C:N ratio
1502               zsurfvar(:,:) = trn(:,:,1,jpphd) * xthetapd
1503#elif defined key_fabm
1504               ! Diatom surface phytoplankton carbon from ERSEM
1505               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_p1c)
1506#else
1507               CALL ctl_stop( ' Trying to run slphydia observation operator', &
1508                  &           ' but no biogeochemical model appears to have been defined' )
1509#endif
1510               llog10 = .TRUE.
1511
1512            CASE('slphynon')
1513#if defined key_hadocc
1514               CALL ctl_stop( ' Trying to run slphynon observation operator', &
1515                  &           ' but HadOCC does not explicitly simulate non-diatoms' )
1516#elif defined key_medusa
1517               ! Non-diatom surface phytoplankton nitrogen from MEDUSA multiplied by C:N ratio
1518               zsurfvar(:,:) = trn(:,:,1,jpphn) * xthetapn
1519#elif defined key_fabm
1520               ! Add all non-diatom surface phytoplankton carbon groups from ERSEM
1521               zsurfvar(:,:) = trn(:,:,1,jp_fabm_m1+jp_fabm_p2c) + &
1522                  &            trn(:,:,1,jp_fabm_m1+jp_fabm_p3c) + trn(:,:,1,jp_fabm_m1+jp_fabm_p4c)
1523#else
1524               CALL ctl_stop( ' Trying to run slphynon observation operator', &
1525                  &           ' but no biogeochemical model appears to have been defined' )
1526#endif
1527               llog10 = .TRUE.
1528
1529            CASE('sspm')
1530#if defined key_spm
1531               zsurfvar(:,:) = 0.0
1532               DO jn = 1, jp_spm
1533                  zsurfvar(:,:) = zsurfvar(:,:) + trn(:,:,1,jn)   ! sum SPM sizes
1534               END DO
1535#else
1536               CALL ctl_stop( ' Trying to run sspm observation operator', &
1537                  &           ' but no spm model appears to have been defined' )
1538#endif
1539
1540            CASE('skd490')
1541#if defined key_hadocc
1542               CALL ctl_stop( ' Trying to run skd490 observation operator', &
1543                  &           ' but HadOCC does not explicitly simulate Kd490' )
1544#elif defined key_medusa
1545               CALL ctl_stop( ' Trying to run skd490 observation operator', &
1546                  &           ' but MEDUSA does not explicitly simulate Kd490' )
1547#elif defined key_fabm
1548               ! light_xEPS diagnostic variable
1549               fabm_3d(:,:,:) = fabm_get_interior_diagnostic_data(model, jp_fabm_xeps)
1550               zsurfvar(:,:) = fabm_3d(:,:,1)
1551#else
1552               CALL ctl_stop( ' Trying to run skd490 observation operator', &
1553                  &           ' but no biogeochemical model appears to have been defined' )
1554#endif
1555
1556            CASE('sfco2')
1557#if defined key_hadocc
1558               zsurfvar(:,:) = HADOCC_FCO2(:,:)    ! fCO2 from HadOCC
1559               IF ( ( MINVAL( HADOCC_FCO2 ) == HADOCC_FILL_FLT ) .AND. &
1560                  & ( MAXVAL( HADOCC_FCO2 ) == HADOCC_FILL_FLT ) ) THEN
1561                  zsurfvar(:,:) = obfillflt
1562                  zsurfmask(:,:) = 0
1563                  CALL ctl_warn( ' HadOCC fCO2 values masked out for observation operator', &
1564                     &           ' as HADOCC_FCO2(:,:) == HADOCC_FILL_FLT' )
1565               ENDIF
1566#elif defined key_medusa && defined key_roam
1567               zsurfvar(:,:) = f2_fco2w(:,:)
1568#elif defined key_fabm
1569               ! First, get pCO2 from FABM
1570               fabm_3d(:,:,:) = fabm_get_interior_diagnostic_data(model, jp_fabm_o3pc)
1571               zsurfvar(:,:) = fabm_3d(:,:,1)
1572               ! Now, convert pCO2 to fCO2, based on SST in K. This follows the standard methodology of:
1573               ! Pierrot et al. (2009), Recommendations for autonomous underway pCO2 measuring systems
1574               ! and data reduction routines, Deep-Sea Research II, 56: 512-522.
1575               ! and
1576               ! Weiss (1974), Carbon dioxide in water and seawater: the solubility of a non-ideal gas,
1577               ! Marine Chemistry, 2: 203-215.
1578               ! In the implementation below, atmospheric pressure has been assumed to be 1 atm and so
1579               ! not explicitly included - atmospheric pressure is not necessarily available so this is
1580               ! the best assumption.
1581               ! Further, the (1-xCO2)^2 term has been neglected. This is common practice
1582               ! (see e.g. Zeebe and Wolf-Gladrow (2001), CO2 in Seawater: Equilibrium, Kinetics, Isotopes)
1583               ! because xCO2 in atm is ~0, and so this term will only affect the result to the 3rd decimal
1584               ! place for typical values, and xCO2 would need to be approximated from pCO2 anyway.
1585               zsurfvar(:,:) = zsurfvar(:,:) * EXP((-1636.75                                                          + &
1586                  &            12.0408      * (tsn(:,:,1,jp_tem)+rt0)                                                 - &
1587                  &            0.0327957    * (tsn(:,:,1,jp_tem)+rt0)*(tsn(:,:,1,jp_tem)+rt0)                         + &
1588                  &            0.0000316528 * (tsn(:,:,1,jp_tem)+rt0)*(tsn(:,:,1,jp_tem)+rt0)*(tsn(:,:,1,jp_tem)+rt0) + &
1589                  &            2.0 * (57.7 - 0.118 * (tsn(:,:,1,jp_tem)+rt0)))                                        / &
1590                  &            (82.0578 * (tsn(:,:,1,jp_tem)+rt0)))
1591#else
1592               CALL ctl_stop( ' Trying to run sfco2 observation operator', &
1593                  &           ' but no biogeochemical model appears to have been defined' )
1594#endif
1595
1596            CASE('spco2')
1597#if defined key_hadocc
1598               zsurfvar(:,:) = HADOCC_PCO2(:,:)    ! pCO2 from HadOCC
1599               IF ( ( MINVAL( HADOCC_PCO2 ) == HADOCC_FILL_FLT ) .AND. &
1600                  & ( MAXVAL( HADOCC_PCO2 ) == HADOCC_FILL_FLT ) ) THEN
1601                  zsurfvar(:,:) = obfillflt
1602                  zsurfmask(:,:) = 0
1603                  CALL ctl_warn( ' HadOCC pCO2 values masked out for observation operator', &
1604                     &           ' as HADOCC_PCO2(:,:) == HADOCC_FILL_FLT' )
1605               ENDIF
1606#elif defined key_medusa && defined key_roam
1607               zsurfvar(:,:) = f2_pco2w(:,:)
1608#elif defined key_fabm
1609               fabm_3d(:,:,:) = fabm_get_interior_diagnostic_data(model, jp_fabm_o3pc)
1610               zsurfvar(:,:) = fabm_3d(:,:,1)
1611#else
1612               CALL ctl_stop( ' Trying to run spco2 observation operator', &
1613                  &           ' but no biogeochemical model appears to have been defined' )
1614#endif
1615
1616            CASE DEFAULT
1617
1618               CALL ctl_stop( 'Unknown surface observation type '//TRIM(cobstypessurf(jtype))//' in dia_obs' )
1619
1620            END SELECT
1621           
1622            IF ( llog10 ) THEN
1623               ! Take the log10 where we can, otherwise exclude
1624               tiny = 1.0e-20
1625               WHERE(zsurfvar(:,:) > tiny .AND. zsurfvar(:,:) /= obfillflt )
1626                  zsurfvar(:,:)  = LOG10(zsurfvar(:,:))
1627               ELSEWHERE
1628                  zsurfvar(:,:)  = obfillflt
1629                  zsurfmask(:,:) = 0
1630               END WHERE
1631            ENDIF
1632
1633            IF ( TRIM(cobstypessurf(jtype)) == 'sla' .AND.                 &
1634                  &  ln_time_mean_sla_bkg ) THEN
1635               !Number of time-steps in meaning period
1636               imeanstp = NINT( ( MeanPeriodHours * 60. * 60. ) / rdt )
1637               CALL obs_surf_opt( surfdataqc(jtype), kstp, jpi, jpj,       &
1638                  &               nit000, idaystp, zsurfvar,               &
1639                  &               zsurfclim, zsurfmask,                    &
1640                  &               n2dintsurf(jtype), llnightav(jtype),     &
1641                  &               ravglamscl(jtype), ravgphiscl(jtype),    &
1642                  &               lfpindegs(jtype), kmeanstp = imeanstp )
1643
1644            ELSE
1645               CALL obs_surf_opt( surfdataqc(jtype), kstp, jpi, jpj,       &
1646                  &               nit000, idaystp, zsurfvar,               &
1647                  &               zsurfclim, zsurfmask,                    &
1648                  &               n2dintsurf(jtype), llnightav(jtype),     &
1649                  &               ravglamscl(jtype), ravgphiscl(jtype),    &
1650                  &               lfpindegs(jtype) )
1651            ENDIF
1652
1653         END DO
1654
1655         CALL wrk_dealloc( jpi, jpj, zsurfvar )
1656         CALL wrk_dealloc( jpi, jpj, zsurfmask )
1657#if defined key_fabm
1658         CALL wrk_dealloc( jpi, jpj, jpk, fabm_3d )
1659#endif
1660
1661      ENDIF
1662
1663   END SUBROUTINE dia_obs
1664
1665   SUBROUTINE dia_obs_wri
1666      !!----------------------------------------------------------------------
1667      !!                    ***  ROUTINE dia_obs_wri  ***
1668      !!         
1669      !! ** Purpose : Call observation diagnostic output routines
1670      !!
1671      !! ** Method  : Call observation diagnostic output routines
1672      !!
1673      !! ** Action  :
1674      !!
1675      !! History :
1676      !!        !  06-03  (K. Mogensen) Original code
1677      !!        !  06-05  (K. Mogensen) Reformatted
1678      !!        !  06-10  (A. Weaver) Cleaning
1679      !!        !  07-03  (K. Mogensen) General handling of profiles
1680      !!        !  08-09  (M. Valdivieso) Velocity component (U,V) profiles
1681      !!        !  15-08  (M. Martin) Combined writing for prof and surf types
1682      !!----------------------------------------------------------------------
1683      !! * Modules used
1684      USE obs_rot_vel          ! Rotation of velocities
1685
1686      IMPLICIT NONE
1687
1688      !! * Local declarations
1689      INTEGER :: jtype                    ! Data set loop variable
1690      INTEGER :: jo, jvar, jk
1691      REAL(wp), DIMENSION(:), ALLOCATABLE :: &
1692         & zu, &
1693         & zv
1694
1695      !-----------------------------------------------------------------------
1696      ! Depending on switches call various observation output routines
1697      !-----------------------------------------------------------------------
1698
1699      IF ( nproftypes > 0 ) THEN
1700
1701         DO jtype = 1, nproftypes
1702
1703            IF ( TRIM(cobstypesprof(jtype)) == 'vel' ) THEN
1704
1705               ! For velocity data, rotate the model velocities to N/S, E/W
1706               ! using the compressed data structure.
1707               ALLOCATE( &
1708                  & zu(profdataqc(jtype)%nvprot(1)), &
1709                  & zv(profdataqc(jtype)%nvprot(2))  &
1710                  & )
1711
1712               CALL obs_rotvel( profdataqc(jtype), nn_2dint_default, zu, zv )
1713
1714               DO jo = 1, profdataqc(jtype)%nprof
1715                  DO jvar = 1, 2
1716                     DO jk = profdataqc(jtype)%npvsta(jo,jvar), profdataqc(jtype)%npvend(jo,jvar)
1717
1718                        IF ( jvar == 1 ) THEN
1719                           profdataqc(jtype)%var(jvar)%vmod(jk) = zu(jk)
1720                        ELSE
1721                           profdataqc(jtype)%var(jvar)%vmod(jk) = zv(jk)
1722                        ENDIF
1723
1724                     END DO
1725                  END DO
1726               END DO
1727
1728               DEALLOCATE( zu )
1729               DEALLOCATE( zv )
1730
1731            END IF
1732
1733            CALL obs_prof_decompress( profdataqc(jtype), &
1734               &                      profdata(jtype), .TRUE., numout )
1735
1736            CALL obs_wri_prof( profdata(jtype) )
1737
1738         END DO
1739
1740      ENDIF
1741
1742      IF ( nsurftypes > 0 ) THEN
1743
1744         DO jtype = 1, nsurftypes
1745
1746            CALL obs_surf_decompress( surfdataqc(jtype), &
1747               &                      surfdata(jtype), .TRUE., numout )
1748
1749            CALL obs_wri_surf( surfdata(jtype) )
1750
1751         END DO
1752
1753      ENDIF
1754
1755   END SUBROUTINE dia_obs_wri
1756
1757   SUBROUTINE dia_obs_dealloc
1758      IMPLICIT NONE
1759      !!----------------------------------------------------------------------
1760      !!                    *** ROUTINE dia_obs_dealloc ***
1761      !!
1762      !!  ** Purpose : To deallocate data to enable the obs_oper online loop.
1763      !!               Specifically: dia_obs_init --> dia_obs --> dia_obs_wri
1764      !!
1765      !!  ** Method : Clean up various arrays left behind by the obs_oper.
1766      !!
1767      !!  ** Action :
1768      !!
1769      !!----------------------------------------------------------------------
1770      ! obs_grid deallocation
1771      CALL obs_grid_deallocate
1772
1773      ! diaobs deallocation
1774      IF ( nproftypes > 0 ) &
1775         &   DEALLOCATE( cobstypesprof, profdata, profdataqc, nvarsprof, nextrprof )
1776
1777      IF ( nsurftypes > 0 ) &
1778         &   DEALLOCATE( cobstypessurf, surfdata, surfdataqc, nvarssurf, nextrsurf, &
1779         &               n2dintsurf, ravglamscl, ravgphiscl, lfpindegs, llnightav )
1780
1781   END SUBROUTINE dia_obs_dealloc
1782
1783   SUBROUTINE ini_date( ddobsini )
1784      !!----------------------------------------------------------------------
1785      !!                    ***  ROUTINE ini_date  ***
1786      !!
1787      !! ** Purpose : Get initial date in double precision YYYYMMDD.HHMMSS format
1788      !!
1789      !! ** Method  : Get initial date in double precision YYYYMMDD.HHMMSS format
1790      !!
1791      !! ** Action  : Get initial date in double precision YYYYMMDD.HHMMSS format
1792      !!
1793      !! History :
1794      !!        !  06-03  (K. Mogensen)  Original code
1795      !!        !  06-05  (K. Mogensen)  Reformatted
1796      !!        !  06-10  (A. Weaver) Cleaning
1797      !!        !  06-10  (G. Smith) Calculates initial date the same as method for final date
1798      !!        !  10-05  (D. Lea) Update to month length calculation for NEMO vn3.2
1799      !!----------------------------------------------------------------------
1800      USE phycst, ONLY : &            ! Physical constants
1801         & rday
1802      USE dom_oce, ONLY : &           ! Ocean space and time domain variables
1803         & rdt
1804
1805      IMPLICIT NONE
1806
1807      !! * Arguments
1808      REAL(dp), INTENT(OUT) :: ddobsini  ! Initial date in YYYYMMDD.HHMMSS
1809
1810      !! * Local declarations
1811      INTEGER :: iyea        ! date - (year, month, day, hour, minute)
1812      INTEGER :: imon
1813      INTEGER :: iday
1814      INTEGER :: ihou
1815      INTEGER :: imin
1816      INTEGER :: imday       ! Number of days in month.
1817      INTEGER, DIMENSION(12) :: &
1818         &       imonth_len  ! Length in days of the months of the current year
1819      REAL(wp) :: zdayfrc    ! Fraction of day
1820
1821      !----------------------------------------------------------------------
1822      ! Initial date initialization (year, month, day, hour, minute)
1823      ! (This assumes that the initial date is for 00z))
1824      !----------------------------------------------------------------------
1825      iyea =   ndate0 / 10000
1826      imon = ( ndate0 - iyea * 10000 ) / 100
1827      iday =   ndate0 - iyea * 10000 - imon * 100
1828      ihou = 0
1829      imin = 0
1830
1831      !----------------------------------------------------------------------
1832      ! Compute number of days + number of hours + min since initial time
1833      !----------------------------------------------------------------------
1834      iday = iday + ( nit000 -1 ) * rdt / rday
1835      zdayfrc = ( nit000 -1 ) * rdt / rday
1836      zdayfrc = zdayfrc - aint(zdayfrc)
1837      ihou = int( zdayfrc * 24 )
1838      imin = int( (zdayfrc * 24 - ihou) * 60 )
1839
1840      !-----------------------------------------------------------------------
1841      ! Convert number of days (iday) into a real date
1842      !----------------------------------------------------------------------
1843
1844      CALL calc_month_len( iyea, imonth_len )
1845
1846      DO WHILE ( iday > imonth_len(imon) )
1847         iday = iday - imonth_len(imon)
1848         imon = imon + 1 
1849         IF ( imon > 12 ) THEN
1850            imon = 1
1851            iyea = iyea + 1
1852            CALL calc_month_len( iyea, imonth_len )  ! update month lengths
1853         ENDIF
1854      END DO
1855
1856      !----------------------------------------------------------------------
1857      ! Convert it into YYYYMMDD.HHMMSS format.
1858      !----------------------------------------------------------------------
1859      ddobsini = iyea * 10000_dp + imon * 100_dp + &
1860         &       iday + ihou * 0.01_dp + imin * 0.0001_dp
1861
1862
1863   END SUBROUTINE ini_date
1864
1865   SUBROUTINE fin_date( ddobsfin )
1866      !!----------------------------------------------------------------------
1867      !!                    ***  ROUTINE fin_date  ***
1868      !!
1869      !! ** Purpose : Get final date in double precision YYYYMMDD.HHMMSS format
1870      !!
1871      !! ** Method  : Get final date in double precision YYYYMMDD.HHMMSS format
1872      !!
1873      !! ** Action  : Get final date in double precision YYYYMMDD.HHMMSS format
1874      !!
1875      !! History :
1876      !!        !  06-03  (K. Mogensen)  Original code
1877      !!        !  06-05  (K. Mogensen)  Reformatted
1878      !!        !  06-10  (A. Weaver) Cleaning
1879      !!        !  10-05  (D. Lea) Update to month length calculation for NEMO vn3.2
1880      !!----------------------------------------------------------------------
1881      USE phycst, ONLY : &            ! Physical constants
1882         & rday
1883      USE dom_oce, ONLY : &           ! Ocean space and time domain variables
1884         & rdt
1885
1886      IMPLICIT NONE
1887
1888      !! * Arguments
1889      REAL(dp), INTENT(OUT) :: ddobsfin ! Final date in YYYYMMDD.HHMMSS
1890
1891      !! * Local declarations
1892      INTEGER :: iyea        ! date - (year, month, day, hour, minute)
1893      INTEGER :: imon
1894      INTEGER :: iday
1895      INTEGER :: ihou
1896      INTEGER :: imin
1897      INTEGER :: imday       ! Number of days in month.
1898      INTEGER, DIMENSION(12) :: &
1899         &       imonth_len  ! Length in days of the months of the current year
1900      REAL(wp) :: zdayfrc    ! Fraction of day
1901
1902      !-----------------------------------------------------------------------
1903      ! Initial date initialization (year, month, day, hour, minute)
1904      ! (This assumes that the initial date is for 00z)
1905      !-----------------------------------------------------------------------
1906      iyea =   ndate0 / 10000
1907      imon = ( ndate0 - iyea * 10000 ) / 100
1908      iday =   ndate0 - iyea * 10000 - imon * 100
1909      ihou = 0
1910      imin = 0
1911
1912      !-----------------------------------------------------------------------
1913      ! Compute number of days + number of hours + min since initial time
1914      !-----------------------------------------------------------------------
1915      iday    = iday +  nitend  * rdt / rday
1916      zdayfrc =  nitend  * rdt / rday
1917      zdayfrc = zdayfrc - AINT( zdayfrc )
1918      ihou    = INT( zdayfrc * 24 )
1919      imin    = INT( ( zdayfrc * 24 - ihou ) * 60 )
1920
1921      !-----------------------------------------------------------------------
1922      ! Convert number of days (iday) into a real date
1923      !----------------------------------------------------------------------
1924
1925      CALL calc_month_len( iyea, imonth_len )
1926
1927      DO WHILE ( iday > imonth_len(imon) )
1928         iday = iday - imonth_len(imon)
1929         imon = imon + 1 
1930         IF ( imon > 12 ) THEN
1931            imon = 1
1932            iyea = iyea + 1
1933            CALL calc_month_len( iyea, imonth_len )  ! update month lengths
1934         ENDIF
1935      END DO
1936
1937      !-----------------------------------------------------------------------
1938      ! Convert it into YYYYMMDD.HHMMSS format
1939      !-----------------------------------------------------------------------
1940      ddobsfin = iyea * 10000_dp + imon * 100_dp    + iday &
1941         &     + ihou * 0.01_dp  + imin * 0.0001_dp
1942
1943    END SUBROUTINE fin_date
1944
1945    SUBROUTINE obs_settypefiles( ntypes, jpmaxnfiles, ifiles, cobstypes, cfiles )
1946
1947       INTEGER, INTENT(IN) :: ntypes      ! Total number of obs types
1948       INTEGER, INTENT(IN) :: jpmaxnfiles ! Maximum number of files allowed for each type
1949       INTEGER, DIMENSION(ntypes), INTENT(OUT) :: &
1950          &                   ifiles      ! Out number of files for each type
1951       CHARACTER(len=8), DIMENSION(ntypes), INTENT(IN) :: &
1952          &                   cobstypes   ! List of obs types
1953       CHARACTER(len=128), DIMENSION(ntypes, jpmaxnfiles), INTENT(IN) :: &
1954          &                   cfiles      ! List of files for all types
1955
1956       !Local variables
1957       INTEGER :: jfile
1958       INTEGER :: jtype
1959
1960       DO jtype = 1, ntypes
1961
1962          ifiles(jtype) = 0
1963          DO jfile = 1, jpmaxnfiles
1964             IF ( trim(cfiles(jtype,jfile)) /= '' ) &
1965                       ifiles(jtype) = ifiles(jtype) + 1
1966          END DO
1967
1968          IF ( ifiles(jtype) == 0 ) THEN
1969               CALL ctl_stop( 'Logical for observation type '//TRIM(cobstypes(jtype))//   &
1970                  &           ' set to true but no files available to read' )
1971          ENDIF
1972
1973          IF(lwp) THEN   
1974             WRITE(numout,*) '             '//cobstypes(jtype)//' input observation file names:'
1975             DO jfile = 1, ifiles(jtype)
1976                WRITE(numout,*) '                '//TRIM(cfiles(jtype,jfile))
1977             END DO
1978          ENDIF
1979
1980       END DO
1981
1982    END SUBROUTINE obs_settypefiles
1983
1984    SUBROUTINE obs_setinterpopts( ntypes, jtype, ctypein,             &
1985               &                  n2dint_default, n2dint_type,        &
1986               &                  ravglamscl_type, ravgphiscl_type,   &
1987               &                  lfp_indegs_type, lavnight_type,     &
1988               &                  n2dint, ravglamscl, ravgphiscl,     &
1989               &                  lfpindegs, lavnight )
1990
1991       INTEGER, INTENT(IN)  :: ntypes             ! Total number of obs types
1992       INTEGER, INTENT(IN)  :: jtype              ! Index of the current type of obs
1993       INTEGER, INTENT(IN)  :: n2dint_default     ! Default option for interpolation type
1994       INTEGER, INTENT(IN)  :: n2dint_type        ! Option for interpolation type
1995       REAL(wp), INTENT(IN) :: &
1996          &                    ravglamscl_type, & !E/W diameter of obs footprint for this type
1997          &                    ravgphiscl_type    !N/S diameter of obs footprint for this type
1998       LOGICAL, INTENT(IN)  :: lfp_indegs_type    !T=> footprint in degrees, F=> in metres
1999       LOGICAL, INTENT(IN)  :: lavnight_type      !T=> obs represent night time average
2000       CHARACTER(len=8), INTENT(IN) :: ctypein 
2001
2002       INTEGER, DIMENSION(ntypes), INTENT(INOUT) :: &
2003          &                    n2dint 
2004       REAL(wp), DIMENSION(ntypes), INTENT(INOUT) :: &
2005          &                    ravglamscl, ravgphiscl
2006       LOGICAL, DIMENSION(ntypes), INTENT(INOUT) :: &
2007          &                    lfpindegs, lavnight
2008
2009       lavnight(jtype) = lavnight_type
2010
2011       IF ( (n2dint_type >= 0) .AND. (n2dint_type <= 6) ) THEN
2012          n2dint(jtype) = n2dint_type
2013       ELSE IF ( n2dint_type == -1 ) THEN
2014          n2dint(jtype) = n2dint_default
2015       ELSE
2016          CALL ctl_stop(' Choice of '//TRIM(ctypein)//' horizontal (2D) interpolation method', &
2017            &                    ' is not available')
2018       ENDIF
2019
2020       ! For averaging observation footprints set options for size of footprint
2021       IF ( (n2dint(jtype) > 4) .AND. (n2dint(jtype) <= 6) ) THEN
2022          IF ( ravglamscl_type > 0._wp ) THEN
2023             ravglamscl(jtype) = ravglamscl_type
2024          ELSE
2025             CALL ctl_stop( 'Incorrect value set for averaging footprint '// &
2026                            'scale (ravglamscl) for observation type '//TRIM(ctypein) )     
2027          ENDIF
2028
2029          IF ( ravgphiscl_type > 0._wp ) THEN
2030             ravgphiscl(jtype) = ravgphiscl_type
2031          ELSE
2032             CALL ctl_stop( 'Incorrect value set for averaging footprint '// &
2033                            'scale (ravgphiscl) for observation type '//TRIM(ctypein) )     
2034          ENDIF
2035
2036          lfpindegs(jtype) = lfp_indegs_type 
2037
2038       ENDIF
2039
2040       ! Write out info
2041       IF(lwp) THEN
2042          IF ( n2dint(jtype) <= 4 ) THEN
2043             WRITE(numout,*) '             '//TRIM(ctypein)// &
2044                &            ' model counterparts will be interpolated horizontally'
2045          ELSE IF ( n2dint(jtype) <= 6 ) THEN
2046             WRITE(numout,*) '             '//TRIM(ctypein)// &
2047                &            ' model counterparts will be averaged horizontally'
2048             WRITE(numout,*) '             '//'    with E/W scale: ',ravglamscl(jtype)
2049             WRITE(numout,*) '             '//'    with N/S scale: ',ravgphiscl(jtype)
2050             IF ( lfpindegs(jtype) ) THEN
2051                 WRITE(numout,*) '             '//'    (in degrees)'
2052             ELSE
2053                 WRITE(numout,*) '             '//'    (in metres)'
2054             ENDIF
2055          ENDIF
2056       ENDIF
2057
2058    END SUBROUTINE obs_setinterpopts
2059
2060END MODULE diaobs
Note: See TracBrowser for help on using the repository browser.