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.
obs_sstbias.F90 in branches/UKMO/dev_r5785_SSS_obsoper/NEMOGCM/NEMO/OPA_SRC/OBS – NEMO

source: branches/UKMO/dev_r5785_SSS_obsoper/NEMOGCM/NEMO/OPA_SRC/OBS/obs_sstbias.F90 @ 7773

Last change on this file since 7773 was 7773, checked in by mattmartin, 7 years ago

Committing updates after doing the following:

  • merging the branch dev_r4650_general_vert_coord_obsoper@7763 into this branch
  • updating it so that the following OBS changes were implemented correctly on top of the simplification changes:
    • generalised vertical coordinate for profile obs. This was done so that is now the default option.
    • sst bias correction implemented with the new simplified obs code.
    • included the biogeochemical obs types int he new simplified obs code.
    • included the changes to exclude obs in the boundary for limited area models
    • included other changes for the efficiency of the obs operator to remove global arrays.
File size: 9.1 KB
Line 
1MODULE obs_sstbias
2   !!======================================================================
3   !!                       ***  MODULE obs_sstbias  ***
4   !! Observation diagnostics: Read the bias for SST data
5   !!======================================================================
6   !!----------------------------------------------------------------------
7   !!   obs_app_sstbias : Driver for reading and applying the SST bias
8   !!----------------------------------------------------------------------
9   !! * Modules used   
10   USE par_kind, ONLY : &       ! Precision variables
11      & wp, &
12      & dp, &
13      & sp
14   USE par_oce, ONLY : &        ! Domain parameters
15      & jpi, &
16      & jpj, &
17      & jpim1
18   USE in_out_manager, ONLY : & ! I/O manager
19      & lwp,    &
20      & numout 
21   USE obs_surf_def             ! Surface observation definitions
22   USE dom_oce, ONLY : &        ! Domain variables
23      & tmask, &
24      & gphit, &
25      & glamt
26   USE obs_inter_h2d
27   USE obs_utils               ! Various observation tools
28   USE obs_inter_sup
29   IMPLICIT NONE
30   !! * Routine accessibility
31   PRIVATE
32   PUBLIC obs_app_sstbias     ! Read the altimeter bias
33CONTAINS
34   SUBROUTINE obs_app_sstbias( sstdata, k2dint, knumtypes, &
35                               cl_bias_files )
36      !!---------------------------------------------------------------------
37      !!
38      !!                   *** ROUTINE obs_app_sstbias ***
39      !!
40      !! ** Purpose : Read SST bias data from files and apply correction to
41      !!             observations
42      !!
43      !! ** Method  :
44      !!
45      !! ** Action  :
46      !!
47      !! References :
48      !!
49      !! History :
50      !!      ! :  2014-08 (J. While) Bias correction code for SST obs,
51      !!      !                      based on obs_rea_altbias
52      !!----------------------------------------------------------------------
53      !! * Modules used
54      USE iom
55      USE netcdf
56
57      !! * Arguments
58      TYPE(obs_surf), INTENT(INOUT) :: &
59         & sstdata            ! SST data
60      INTEGER, INTENT(IN) :: k2dint
61      INTEGER, INTENT(IN) :: &
62         & knumtypes          ! Number of bias types to read in
63      CHARACTER(LEN=128), DIMENSION(knumtypes), INTENT(IN) :: &
64         & cl_bias_files      ! List of files to read
65
66      !! * Local declarations
67      INTEGER :: jslano       ! Data set loop variable
68      INTEGER :: jobs         ! Obs loop variable
69      INTEGER :: jpisstbias   ! Number of grid point in latitude for the bias
70      INTEGER :: jpjsstbias   ! Number of grid point in longitude for the bias
71      INTEGER :: iico         ! Grid point indices
72      INTEGER :: ijco
73      INTEGER :: jt
74      INTEGER :: i_nx_id      ! Index to read the NetCDF file
75      INTEGER :: i_ny_id      !
76      INTEGER :: i_file_id    !
77      INTEGER :: i_var_id
78      INTEGER, DIMENSION(knumtypes) :: &
79         & ibiastypes         ! Array of the bias types in each file
80      REAL(wp), DIMENSION(jpi,jpj,knumtypes) :: & 
81         & z_sstbias          ! Array to store the SST bias values
82      REAL(wp), DIMENSION(jpi,jpj) :: & 
83         & z_sstbias_2d       ! Array to store the SST bias values   
84      REAL(wp), DIMENSION(1) :: &
85         & zext, &
86         & zobsmask
87      REAL(wp), DIMENSION(2,2,1) :: &
88         & zweig
89      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: &
90         & zmask, &
91         & zglam, &
92         & zgphi
93      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: &
94         & zmask_tmp, &
95         & zglam_tmp, &
96         & zgphi_tmp   
97      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::  zbias   
98      REAL(wp) :: zlam
99      REAL(wp) :: zphi
100      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: &
101         & igrdi, &
102         & igrdj
103      INTEGER, DIMENSION(:,:,:), ALLOCATABLE :: &
104         & igrdi_tmp, &
105         & igrdj_tmp   
106      INTEGER :: numsstbias
107      INTEGER(KIND=NF90_INT) :: ifile_source
108     
109      INTEGER :: incfile
110      INTEGER :: jtype
111      INTEGER :: iret 
112      INTEGER :: inumtype
113
114      IF ( lwp ) THEN
115         WRITE(numout,*) 
116         WRITE(numout,*) 'obs_app_sstbias : '
117         WRITE(numout,*) '----------------- '
118         WRITE(numout,*) 'Read SST bias '
119      ENDIF
120
121      ! Open and read the SST bias files for each bias type
122      z_sstbias(:,:,:) = 0.0_wp
123
124      DO jtype = 1, knumtypes
125     
126         numsstbias=0
127
128         IF ( lwp ) WRITE(numout,*) 'Opening ',cl_bias_files(jtype)
129         CALL iom_open( cl_bias_files(jtype), numsstbias, ldstop=.FALSE. )
130
131         IF (numsstbias .GT. 0) THEN
132     
133            !Read the bias type from the file
134            !No IOM get attribute command at time of writing,
135            !so have to use NETCDF
136            !routines directly - should be upgraded in the future
137            iret=NF90_OPEN(TRIM(cl_bias_files(jtype)), NF90_NOWRITE, incfile)
138            iret=NF90_GET_ATT( incfile, NF90_GLOBAL, "SST_source", &
139                              ifile_source )
140            ibiastypes(jtype) = ifile_source                 
141            iret=NF90_CLOSE(incfile)       
142           
143            IF ( iret /= 0  ) THEN
144               CALL ctl_stop( 'obs_app_sstbias : Cannot read bias type from file '// &
145                  &           TRIM( cl_bias_files(jtype) ) )
146            ENDIF
147
148            ! Get the SST bias data
149            CALL iom_get( numsstbias, jpdom_data, 'tn', z_sstbias_2d(:,:), 1 )
150            z_sstbias(:,:,jtype) = z_sstbias_2d(:,:)       
151            ! Close the file
152            CALL iom_close(numsstbias)
153     
154         ELSE
155            CALL ctl_stop('obs_read_sstbias: File '// & 
156               &          TRIM( cl_bias_files(jtype) )//' Not found')
157         ENDIF
158
159      END DO
160           
161      ! Interpolate the bias from the model grid to the observation points
162      ALLOCATE( &
163         & igrdi(2,2,sstdata%nsurf), &
164         & igrdj(2,2,sstdata%nsurf), &
165         & zglam(2,2,sstdata%nsurf), &
166         & zgphi(2,2,sstdata%nsurf), &
167         & zmask(2,2,sstdata%nsurf)  )
168       
169      DO jobs = 1, sstdata%nsurf 
170         igrdi(1,1,jobs) = sstdata%mi(jobs)-1
171         igrdj(1,1,jobs) = sstdata%mj(jobs)-1
172         igrdi(1,2,jobs) = sstdata%mi(jobs)-1
173         igrdj(1,2,jobs) = sstdata%mj(jobs)
174         igrdi(2,1,jobs) = sstdata%mi(jobs)
175         igrdj(2,1,jobs) = sstdata%mj(jobs)-1
176         igrdi(2,2,jobs) = sstdata%mi(jobs)
177         igrdj(2,2,jobs) = sstdata%mj(jobs)
178      END DO
179
180      CALL obs_int_comm_2d( 2, 2, sstdata%nsurf, &
181         &                  igrdi, igrdj, glamt, zglam )
182      CALL obs_int_comm_2d( 2, 2, sstdata%nsurf, &
183         &                  igrdi, igrdj, gphit, zgphi )
184      CALL obs_int_comm_2d( 2, 2, sstdata%nsurf, &
185         &                  igrdi, igrdj, tmask(:,:,1), zmask )
186
187      DO jtype = 1, knumtypes
188         
189         !Find the number observations of type
190         !and alllocate tempory arrays
191         inumtype = COUNT( sstdata%ntyp(:) == ibiastypes(jtype) )
192
193         ALLOCATE( &
194            & igrdi_tmp(2,2,inumtype), &
195            & igrdj_tmp(2,2,inumtype), &
196            & zglam_tmp(2,2,inumtype), &
197            & zgphi_tmp(2,2,inumtype), &
198            & zmask_tmp(2,2,inumtype), &
199            & zbias( 2,2,inumtype ) )
200
201         jt=1
202         DO jobs = 1, sstdata%nsurf 
203            IF ( sstdata%ntyp(jobs) == ibiastypes(jtype) ) THEN
204
205               igrdi_tmp(:,:,jt) = igrdi(:,:,jobs) 
206               igrdj_tmp(:,:,jt) = igrdj(:,:,jobs)
207               zglam_tmp(:,:,jt) = zglam(:,:,jobs)
208               zgphi_tmp(:,:,jt) = zgphi(:,:,jobs)
209               zgphi_tmp(:,:,jt) = zgphi(:,:,jobs)
210               zmask_tmp(:,:,jt) = zmask(:,:,jobs)
211
212               jt = jt +1
213
214            ENDIF
215         END DO
216                         
217         CALL obs_int_comm_2d( 2, 2, inumtype, &
218               &           igrdi_tmp(:,:,:), igrdj_tmp(:,:,:), &
219               &           z_sstbias(:,:,jtype), zbias(:,:,:) )
220
221         jt=1
222         DO jobs = 1, sstdata%nsurf
223            IF ( sstdata%ntyp(jobs) == ibiastypes(jtype) ) THEN
224
225               zlam = sstdata%rlam(jobs)
226               zphi = sstdata%rphi(jobs)
227               iico = sstdata%mi(jobs)
228               ijco = sstdata%mj(jobs)         
229
230               CALL obs_int_h2d_init( 1, 1, k2dint, zlam, zphi,         &
231                  &                   zglam_tmp(:,:,jt), &
232                  &                   zgphi_tmp(:,:,jt), &
233                  &                   zmask_tmp(:,:,jt), zweig, zobsmask )
234
235               CALL obs_int_h2d( 1, 1, zweig, zbias(:,:,jt),  zext )
236
237               ! adjust sst with bias field
238               sstdata%robs(jobs,1) = &
239                  &    sstdata%robs(jobs,1) - zext(1)
240
241               jt=jt+1
242
243            ENDIF
244         END DO 
245               
246         !Deallocate arrays
247         DEALLOCATE( &
248         & igrdi_tmp, &
249         & igrdj_tmp, &
250         & zglam_tmp, &
251         & zgphi_tmp, &
252         & zmask_tmp, &
253         & zbias )     
254     
255      END DO !jtype
256
257      DEALLOCATE( &
258         & igrdi, &
259         & igrdj, &
260         & zglam, &
261         & zgphi, &
262         & zmask )
263
264      IF(lwp) THEN
265         WRITE(numout,*) " "
266         WRITE(numout,*) "SST bias correction applied successfully"
267         WRITE(numout,*) "Obs types: ",ibiastypes(:), &
268                              " have all been bias corrected\n"
269      ENDIF
270   END SUBROUTINE obs_app_sstbias
271 
272END MODULE obs_sstbias
Note: See TracBrowser for help on using the repository browser.