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.
fldread.F90 in NEMO/branches/2019/dev_r11351_fldread_with_XIOS/src/OCE/SBC – NEMO

source: NEMO/branches/2019/dev_r11351_fldread_with_XIOS/src/OCE/SBC/fldread.F90 @ 11482

Last change on this file since 11482 was 11482, checked in by andmirek, 5 years ago

Ticket #2195 read initial conditions with XIOS

  • Property svn:keywords set to Id
File size: 110.5 KB
Line 
1MODULE fldread
2   !!======================================================================
3   !!                       ***  MODULE  fldread  ***
4   !! Ocean forcing:  read input field for surface boundary condition
5   !!=====================================================================
6   !! History :  2.0  !  2006-06  (S. Masson, G. Madec)  Original code
7   !!            3.0  !  2008-05  (S. Alderson)  Modified for Interpolation in memory from input grid to model grid
8   !!            3.4  !  2013-10  (D. Delrosso, P. Oddo)  suppression of land point prior to interpolation
9   !!                 !  12-2015  (J. Harle) Adding BDY on-the-fly interpolation
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   fld_read      : read input fields used for the computation of the surface boundary condition
14   !!   fld_init      : initialization of field read
15   !!   fld_rec       : determined the record(s) to be read
16   !!   fld_get       : read the data
17   !!   fld_map       : read global data from file and map onto local data using a general mapping (use for open boundaries)
18   !!   fld_rot       : rotate the vector fields onto the local grid direction
19   !!   fld_clopn     : update the data file name and close/open the files
20   !!   fld_fill      : fill the data structure with the associated information read in namelist
21   !!   wgt_list      : manage the weights used for interpolation
22   !!   wgt_print     : print the list of known weights
23   !!   fld_weight    : create a WGT structure and fill in data from file, restructuring as required
24   !!   apply_seaoverland : fill land with ocean values
25   !!   seaoverland   : create shifted matrices for seaoverland application
26   !!   fld_interp    : apply weights to input gridded data to create data on model grid
27   !!   ksec_week     : function returning the first 3 letters of the first day of the weekly file
28   !!----------------------------------------------------------------------
29   USE oce            ! ocean dynamics and tracers
30   USE dom_oce        ! ocean space and time domain
31   USE phycst         ! physical constant
32   USE sbc_oce        ! surface boundary conditions : fields
33   USE geo2ocean      ! for vector rotation on to model grid
34   !
35   USE in_out_manager ! I/O manager
36   USE iom            ! I/O manager library
37   USE ioipsl  , ONLY : ymds2ju, ju2ymds   ! for calendar
38   USE lib_mpp        ! MPP library
39   USE lbclnk         ! ocean lateral boundary conditions (C1D case)
40   
41   IMPLICIT NONE
42   PRIVATE   
43 
44   PUBLIC   fld_map    ! routine called by tides_init
45   PUBLIC   fld_read, fld_fill   ! called by sbc... modules
46   PUBLIC   fld_clopn
47
48   TYPE, PUBLIC ::   FLD_N      !: Namelist field informations
49      CHARACTER(len = 256) ::   clname      ! generic name of the NetCDF flux file
50      REAL(wp)             ::   nfreqh      ! frequency of each flux file
51      CHARACTER(len = 34)  ::   clvar       ! generic name of the variable in the NetCDF flux file
52      LOGICAL              ::   ln_tint     ! time interpolation or not (T/F)
53      LOGICAL              ::   ln_clim     ! climatology or not (T/F)
54      CHARACTER(len = 8)   ::   cltype      ! type of data file 'daily', 'monthly' or yearly'
55      CHARACTER(len = 256) ::   wname       ! generic name of a NetCDF weights file to be used, blank if not
56      CHARACTER(len = 34)  ::   vcomp       ! symbolic component name if a vector that needs rotation
57      !                                     ! a string starting with "U" or "V" for each component   
58      !                                     ! chars 2 onwards identify which components go together 
59      CHARACTER(len = 34)  ::   lname       ! generic name of a NetCDF land/sea mask file to be used, blank if not
60      !                                     ! 0=sea 1=land
61   END TYPE FLD_N
62
63   TYPE, PUBLIC ::   FLD        !: Input field related variables
64      CHARACTER(len = 256)            ::   clrootname   ! generic name of the NetCDF file
65      CHARACTER(len = 256)            ::   clname       ! current name of the NetCDF file
66      REAL(wp)                        ::   nfreqh       ! frequency of each flux file
67      CHARACTER(len = 34)             ::   clvar        ! generic name of the variable in the NetCDF flux file
68      LOGICAL                         ::   ln_tint      ! time interpolation or not (T/F)
69      LOGICAL                         ::   ln_clim      ! climatology or not (T/F)
70      CHARACTER(len = 8)              ::   cltype       ! type of data file 'daily', 'monthly' or yearly'
71      INTEGER                         ::   num          ! iom id of the jpfld files to be read
72      INTEGER , DIMENSION(2)          ::   nrec_b       ! before record (1: index, 2: second since Jan. 1st 00h of nit000 year)
73      INTEGER , DIMENSION(2)          ::   nrec_a       ! after  record (1: index, 2: second since Jan. 1st 00h of nit000 year)
74      REAL(wp) , ALLOCATABLE, DIMENSION(:,:,:  ) ::   fnow   ! input fields interpolated to now time step
75      REAL(wp) , ALLOCATABLE, DIMENSION(:,:,:,:) ::   fdta   ! 2 consecutive record of input fields
76      CHARACTER(len = 256)            ::   wgtname      ! current name of the NetCDF weight file acting as a key
77      !                                                 ! into the WGTLIST structure
78      CHARACTER(len = 34)             ::   vcomp        ! symbolic name for a vector component that needs rotation
79      LOGICAL, DIMENSION(2)           ::   rotn         ! flag to indicate whether before/after field has been rotated
80      INTEGER                         ::   nreclast     ! last record to be read in the current file
81      CHARACTER(len = 256)            ::   lsmname      ! current name of the NetCDF mask file acting as a key
82      INTEGER                         ::   igrd         ! grid type for bdy data
83      INTEGER                         ::   ibdy         ! bdy set id number
84   END TYPE FLD
85
86   TYPE, PUBLIC ::   MAP_POINTER      !: Map from input data file to local domain
87      INTEGER, POINTER, DIMENSION(:)  ::  ptr           ! Array of integer pointers to 1D arrays
88      LOGICAL                         ::  ll_unstruc    ! Unstructured (T) or structured (F) boundary data file
89   END TYPE MAP_POINTER
90
91!$AGRIF_DO_NOT_TREAT
92
93   !! keep list of all weights variables so they're only read in once
94   !! need to add AGRIF directives not to process this structure
95   !! also need to force wgtname to include AGRIF nest number
96   TYPE         ::   WGT        !: Input weights related variables
97      CHARACTER(len = 256)                    ::   wgtname      ! current name of the NetCDF weight file
98      INTEGER , DIMENSION(2)                  ::   ddims        ! shape of input grid
99      INTEGER , DIMENSION(2)                  ::   botleft      ! top left corner of box in input grid containing
100      !                                                         ! current processor grid
101      INTEGER , DIMENSION(2)                  ::   topright     ! top right corner of box
102      INTEGER                                 ::   jpiwgt       ! width of box on input grid
103      INTEGER                                 ::   jpjwgt       ! height of box on input grid
104      INTEGER                                 ::   numwgt       ! number of weights (4=bilinear, 16=bicubic)
105      INTEGER                                 ::   nestid       ! for agrif, keep track of nest we're in
106      INTEGER                                 ::   overlap      ! =0 when cyclic grid has no overlapping EW columns
107      !                                                         ! =>1 when they have one or more overlapping columns     
108      !                                                         ! =-1 not cyclic
109      LOGICAL                                 ::   cyclic       ! east-west cyclic or not
110      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpi     ! array of source integers
111      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpj     ! array of source integers
112      REAL(wp), DIMENSION(:,:,:), POINTER     ::   data_wgt     ! array of weights on model grid
113      REAL(wp), DIMENSION(:,:,:), POINTER     ::   fly_dta      ! array of values on input grid
114      REAL(wp), DIMENSION(:,:,:), POINTER     ::   col          ! temporary array for reading in columns
115   END TYPE WGT
116
117   INTEGER,     PARAMETER             ::   tot_wgts = 20
118   TYPE( WGT ), DIMENSION(tot_wgts)   ::   ref_wgts     ! array of wgts
119   INTEGER                            ::   nxt_wgt = 1  ! point to next available space in ref_wgts array
120   REAL(wp), PARAMETER                ::   undeff_lsm = -999.00_wp
121
122!$AGRIF_END_DO_NOT_TREAT
123
124   !!----------------------------------------------------------------------
125   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
126   !! $Id$
127   !! Software governed by the CeCILL license (see ./LICENSE)
128   !!----------------------------------------------------------------------
129CONTAINS
130
131   SUBROUTINE fld_read( kt, kn_fsbc, sd, map, kit, kt_offset, jpk_bdy, fvl, ldxios )
132      !!---------------------------------------------------------------------
133      !!                    ***  ROUTINE fld_read  ***
134      !!                   
135      !! ** Purpose :   provide at each time step the surface ocean fluxes
136      !!                (momentum, heat, freshwater and runoff)
137      !!
138      !! ** Method  :   READ each input fields in NetCDF files using IOM
139      !!      and intepolate it to the model time-step.
140      !!         Several assumptions are made on the input file:
141      !!      blahblahblah....
142      !!----------------------------------------------------------------------
143      INTEGER  , INTENT(in   )               ::   kt        ! ocean time step
144      INTEGER  , INTENT(in   )               ::   kn_fsbc   ! sbc computation period (in time step)
145      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
146      TYPE(MAP_POINTER),INTENT(in), OPTIONAL, DIMENSION(:) ::   map   ! global-to-local mapping indices
147      INTEGER  , INTENT(in   ), OPTIONAL     ::   kit       ! subcycle timestep for timesplitting option
148      INTEGER  , INTENT(in   ), OPTIONAL     ::   kt_offset ! provide fields at time other than "now"
149      !                                                     !   kt_offset = -1 => fields at "before" time level
150      !                                                     !   kt_offset = +1 => fields at "after"  time level
151      !                                                     !   etc.
152      INTEGER  , INTENT(in   ), OPTIONAL     ::   jpk_bdy   ! number of vertical levels in the BDY data
153      LOGICAL  , INTENT(in   ), OPTIONAL     ::   fvl   ! number of vertical levels in the BDY data
154      LOGICAL  , INTENT(in   ), OPTIONAL     ::   ldxios    ! use xios to read data
155      !!
156      INTEGER  ::   itmp         ! local variable
157      INTEGER  ::   imf          ! size of the structure sd
158      INTEGER  ::   jf           ! dummy indices
159      INTEGER  ::   isecend      ! number of second since Jan. 1st 00h of nit000 year at nitend
160      INTEGER  ::   isecsbc      ! number of seconds between Jan. 1st 00h of nit000 year and the middle of sbc time step
161      INTEGER  ::   it_offset    ! local time offset variable
162      LOGICAL  ::   llnxtyr      ! open next year  file?
163      LOGICAL  ::   llnxtmth     ! open next month file?
164      LOGICAL  ::   llstop       ! stop is the file does not exist
165      LOGICAL  ::   ll_firstcall ! true if this is the first call to fld_read for this set of fields
166      REAL(wp) ::   ztinta       ! ratio applied to after  records when doing time interpolation
167      REAL(wp) ::   ztintb       ! ratio applied to before records when doing time interpolation
168      CHARACTER(LEN=1000) ::   clfmt  ! write format
169      TYPE(MAP_POINTER)   ::   imap   ! global-to-local mapping indices
170      LOGICAL  ::   luxios       ! temporary variable to hold ldxios
171      !!---------------------------------------------------------------------
172
173      luxios = .FALSE.
174      IF( PRESENT(ldxios) ) luxios = ldxios
175
176      ll_firstcall = kt == nit000
177      IF( PRESENT(kit) )   ll_firstcall = ll_firstcall .and. kit == 1
178
179      IF ( nn_components == jp_iam_sas ) THEN   ;   it_offset = nn_fsbc
180      ELSE                                      ;   it_offset = 0
181      ENDIF
182      IF( PRESENT(kt_offset) )   it_offset = kt_offset
183
184      imap%ptr => NULL()
185
186      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
187      IF( present(kit) ) THEN   ! ignore kn_fsbc in this case
188         isecsbc = nsec_year + nsec1jan000 + (kit+it_offset)*NINT( rdt/REAL(nn_baro,wp) )
189      ELSE                      ! middle of sbc time step
190         isecsbc = nsec_year + nsec1jan000 + NINT(0.5 * REAL(kn_fsbc - 1,wp) * rdt) + it_offset * NINT(rdt)
191      ENDIF
192      imf = SIZE( sd )
193      !
194      IF( ll_firstcall ) THEN                      ! initialization
195         DO jf = 1, imf 
196            IF( PRESENT(map) ) imap = map(jf)
197               IF( PRESENT(jpk_bdy) ) THEN
198                  CALL fld_init( kn_fsbc, sd(jf), imap, jpk_bdy, fvl, ldxios = luxios)  ! read each before field (put them in after as they will be swapped)
199               ELSE
200                  CALL fld_init( kn_fsbc, sd(jf), imap, ldxios = luxios)  ! read each before field (put them in after as they will be swapped)
201               ENDIF
202         END DO
203         IF( lwp ) CALL wgt_print()                ! control print
204      ENDIF
205      !                                            ! ====================================== !
206      IF( MOD( kt-1, kn_fsbc ) == 0 ) THEN         ! update field at each kn_fsbc time-step !
207         !                                         ! ====================================== !
208         !
209         DO jf = 1, imf                            ! ---   loop over field   --- !
210           
211            IF( isecsbc > sd(jf)%nrec_a(2) .OR. ll_firstcall ) THEN    ! read/update the after data?
212
213               IF( PRESENT(map) )   imap = map(jf)   ! temporary definition of map
214
215               sd(jf)%nrec_b(:) = sd(jf)%nrec_a(:)                                  ! swap before record informations
216               sd(jf)%rotn(1) = sd(jf)%rotn(2)                                      ! swap before rotate informations
217               IF( sd(jf)%ln_tint )   sd(jf)%fdta(:,:,:,1) = sd(jf)%fdta(:,:,:,2)   ! swap before record field
218
219               CALL fld_rec( kn_fsbc, sd(jf), kt_offset = it_offset, kit = kit)    ! update after record informations
220
221               ! if kn_fsbc*rdt is larger than nfreqh (which is kind of odd),
222               ! it is possible that the before value is no more the good one... we have to re-read it
223               ! if before is not the last record of the file currently opened and after is the first record to be read
224               ! in a new file which means after = 1 (the file to be opened corresponds to the current time)
225               ! or after = nreclast + 1 (the file to be opened corresponds to a future time step)
226               IF( .NOT. ll_firstcall .AND. sd(jf)%ln_tint .AND. sd(jf)%nrec_b(1) /= sd(jf)%nreclast &
227                  &                   .AND. MOD( sd(jf)%nrec_a(1), sd(jf)%nreclast ) == 1 ) THEN
228                  itmp = sd(jf)%nrec_a(1)                       ! temporary storage
229                  sd(jf)%nrec_a(1) = sd(jf)%nreclast            ! read the last record of the file currently opened
230                  CALL fld_get( sd(jf), imap, ldxios = luxios)  ! read after data
231                  sd(jf)%fdta(:,:,:,1) = sd(jf)%fdta(:,:,:,2)   ! re-swap before record field
232                  sd(jf)%nrec_b(1) = sd(jf)%nrec_a(1)           ! update before record informations
233                  sd(jf)%nrec_b(2) = sd(jf)%nrec_a(2) - NINT( sd(jf)%nfreqh * 3600 )  ! assume freq to be in hours in this case
234                  sd(jf)%rotn(1)   = sd(jf)%rotn(2)             ! update before rotate informations
235                  sd(jf)%nrec_a(1) = itmp                       ! move back to after record
236               ENDIF
237
238               CALL fld_clopn( sd(jf) )   ! Do we need to open a new year/month/week/day file?
239               
240               IF( sd(jf)%ln_tint ) THEN
241                 
242                  ! if kn_fsbc*rdt is larger than nfreqh (which is kind of odd),
243                  ! it is possible that the before value is no more the good one... we have to re-read it
244                  ! if before record is not just just before the after record...
245                  IF( .NOT. ll_firstcall .AND. MOD( sd(jf)%nrec_a(1), sd(jf)%nreclast ) /= 1 &
246                     &                   .AND. sd(jf)%nrec_b(1) /= sd(jf)%nrec_a(1) - 1 ) THEN   
247                     sd(jf)%nrec_a(1) = sd(jf)%nrec_a(1) - 1       ! move back to before record
248                     CALL fld_get( sd(jf), imap, ldxios = luxios)  ! read after data
249                     sd(jf)%fdta(:,:,:,1) = sd(jf)%fdta(:,:,:,2)   ! re-swap before record field
250                     sd(jf)%nrec_b(1) = sd(jf)%nrec_a(1)           ! update before record informations
251                     sd(jf)%nrec_b(2) = sd(jf)%nrec_a(2) - NINT( sd(jf)%nfreqh * 3600 )  ! assume freq to be in hours in this case
252                     sd(jf)%rotn(1)   = sd(jf)%rotn(2)             ! update before rotate informations
253                     sd(jf)%nrec_a(1) = sd(jf)%nrec_a(1) + 1       ! move back to after record
254                  ENDIF
255               ENDIF ! temporal interpolation?
256
257               ! do we have to change the year/month/week/day of the forcing field??
258               ! if we do time interpolation we will need to open next year/month/week/day file before the end of the current
259               ! one. If so, we are still before the end of the year/month/week/day when calling fld_rec so sd(jf)%nrec_a(1)
260               ! will be larger than the record number that should be read for current year/month/week/day
261               ! do we need next file data?
262               ! This applies to both cases with or without time interpolation
263               IF( sd(jf)%nrec_a(1) > sd(jf)%nreclast ) THEN
264                 
265                  sd(jf)%nrec_a(1) = sd(jf)%nrec_a(1) - sd(jf)%nreclast   !
266                 
267                  IF( .NOT. ( sd(jf)%ln_clim .AND. sd(jf)%cltype == 'yearly' ) ) THEN   ! close/open the current/new file
268                     
269                     llnxtmth = sd(jf)%cltype == 'monthly' .OR. nday == nmonth_len(nmonth)      ! open next month file?
270                     llnxtyr  = sd(jf)%cltype == 'yearly'  .OR. (nmonth == 12 .AND. llnxtmth)   ! open next year  file?
271
272                     ! if the run finishes at the end of the current year/month/week/day, we will allow next
273                     ! year/month/week/day file to be not present. If the run continue further than the current
274                     ! year/month/week/day, next year/month/week/day file must exist
275                     isecend = nsec_year + nsec1jan000 + (nitend - kt) * NINT(rdt)   ! second at the end of the run
276                     llstop = isecend > sd(jf)%nrec_a(2)                                   ! read more than 1 record of next year
277                     ! we suppose that the date of next file is next day (should be ok even for weekly files...)
278                     CALL fld_clopn( sd(jf), nyear  + COUNT((/llnxtyr /))                                           ,         &
279                        &                    nmonth + COUNT((/llnxtmth/)) - 12                 * COUNT((/llnxtyr /)),         &
280                        &                    nday   + 1                   - nmonth_len(nmonth) * COUNT((/llnxtmth/)), llstop )
281
282                     IF( sd(jf)%num <= 0 .AND. .NOT. llstop ) THEN    ! next year file does not exist
283                        CALL ctl_warn('next year/month/week/day file: '//TRIM(sd(jf)%clname)//     &
284                           &     ' not present -> back to current year/month/day')
285                        CALL fld_clopn( sd(jf) )       ! back to the current year/month/day
286                        sd(jf)%nrec_a(1) = sd(jf)%nreclast     ! force to read the last record in the current year file
287                     ENDIF
288                     
289                  ENDIF
290               ENDIF   ! open need next file?
291                 
292               ! read after data
293               IF( PRESENT(jpk_bdy) ) THEN
294                  CALL fld_get( sd(jf), imap, jpk_bdy, fvl, ldxios = luxios)
295               ELSE
296                  CALL fld_get( sd(jf), imap, ldxios = luxios )
297               ENDIF
298            ENDIF   ! read new data?
299         END DO                                    ! --- end loop over field --- !
300
301         CALL fld_rot( kt, sd )                    ! rotate vector before/now/after fields if needed
302
303         DO jf = 1, imf                            ! ---   loop over field   --- !
304            !
305            IF( sd(jf)%ln_tint ) THEN              ! temporal interpolation
306               IF(lwp .AND. kt - nit000 <= 100 ) THEN
307                  clfmt = "('   fld_read: var ', a, ' kt = ', i8, ' (', f9.4,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
308                     &    "', records b/a: ', i6.4, '/', i6.4, ' (days ', f9.4,'/', f9.4, ')')"
309                  WRITE(numout, clfmt)  TRIM( sd(jf)%clvar ), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,   &           
310                     & sd(jf)%nrec_b(1), sd(jf)%nrec_a(1), REAL(sd(jf)%nrec_b(2),wp)/rday, REAL(sd(jf)%nrec_a(2),wp)/rday
311                  WRITE(numout, *) '      it_offset is : ',it_offset
312               ENDIF
313               ! temporal interpolation weights
314               ztinta =  REAL( isecsbc - sd(jf)%nrec_b(2), wp ) / REAL( sd(jf)%nrec_a(2) - sd(jf)%nrec_b(2), wp )
315               ztintb =  1. - ztinta
316               sd(jf)%fnow(:,:,:) = ztintb * sd(jf)%fdta(:,:,:,1) + ztinta * sd(jf)%fdta(:,:,:,2)
317            ELSE   ! nothing to do...
318               IF(lwp .AND. kt - nit000 <= 100 ) THEN
319                  clfmt = "('   fld_read: var ', a, ' kt = ', i8,' (', f9.4,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
320                     &    "', record: ', i6.4, ' (days ', f9.4, ' <-> ', f9.4, ')')"
321                  WRITE(numout, clfmt) TRIM(sd(jf)%clvar), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,    &
322                     &                 sd(jf)%nrec_a(1), REAL(sd(jf)%nrec_b(2),wp)/rday, REAL(sd(jf)%nrec_a(2),wp)/rday
323               ENDIF
324            ENDIF
325            !
326            IF( kt == nitend - kn_fsbc + 1 )   CALL iom_close( sd(jf)%num )   ! Close the input files
327
328         END DO                                    ! --- end loop over field --- !
329         !
330      ENDIF
331      !
332   END SUBROUTINE fld_read
333
334
335   SUBROUTINE fld_init( kn_fsbc, sdjf, map , jpk_bdy, fvl, ldxios)
336      !!---------------------------------------------------------------------
337      !!                    ***  ROUTINE fld_init  ***
338      !!
339      !! ** Purpose :  - first call to fld_rec to define before values
340      !!               - if time interpolation, read before data
341      !!----------------------------------------------------------------------
342      INTEGER  , INTENT(in   ) ::   kn_fsbc      ! sbc computation period (in time step)
343      TYPE(FLD), INTENT(inout) ::   sdjf         ! input field related variables
344      TYPE(MAP_POINTER),INTENT(in) ::   map      ! global-to-local mapping indices
345      INTEGER  , INTENT(in), OPTIONAL :: jpk_bdy ! number of vertical levels in the BDY data
346      LOGICAL  , INTENT(in), OPTIONAL :: fvl     ! number of vertical levels in the BDY data
347      LOGICAL  ,             OPTIONAL :: ldxios  ! use xios for I/O
348      !!
349      LOGICAL :: llprevyr              ! are we reading previous year  file?
350      LOGICAL :: llprevmth             ! are we reading previous month file?
351      LOGICAL :: llprevweek            ! are we reading previous week  file?
352      LOGICAL :: llprevday             ! are we reading previous day   file?
353      LOGICAL :: llprev                ! llprevyr .OR. llprevmth .OR. llprevweek .OR. llprevday
354      INTEGER :: idvar                 ! variable id
355      INTEGER :: inrec                 ! number of record existing for this variable
356      INTEGER :: iyear, imonth, iday   ! first day of the current file in yyyy mm dd
357      INTEGER :: isec_week             ! number of seconds since start of the weekly file
358      CHARACTER(LEN=1000) ::   clfmt   ! write format
359      LOGICAL :: luxios                ! local variable to keep ldxios info
360      !!---------------------------------------------------------------------
361      llprevyr   = .FALSE.
362      llprevmth  = .FALSE.
363      llprevweek = .FALSE.
364      llprevday  = .FALSE.
365      isec_week  = 0
366      luxios = .FALSE.
367      IF( PRESENT(ldxios)) luxios = ldxios
368
369      !
370      ! define record informations
371      CALL fld_rec( kn_fsbc, sdjf, ldbefore = .TRUE. )  ! return before values in sdjf%nrec_a (as we will swap it later)
372      !
373      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
374      !
375      IF( sdjf%ln_tint ) THEN ! we need to read the previous record and we will put it in the current record structure
376         !
377         IF( sdjf%nrec_a(1) == 0  ) THEN   ! we redefine record sdjf%nrec_a(1) with the last record of previous year file
378            IF(luxios) THEN
379               IF(lwp) write(numout, *) 'Can not use XIOS to read data. Switching to old functionality' 
380               luxios = .FALSE.               ! we can not use XIOS in this case
381            ENDIF
382            IF    ( sdjf%nfreqh == -12 ) THEN   ! yearly mean
383               IF( sdjf%cltype == 'yearly' ) THEN             ! yearly file
384                  sdjf%nrec_a(1) = 1                                                       ! force to read the unique record
385                  llprevyr  = .NOT. sdjf%ln_clim                                           ! use previous year  file?
386               ELSE
387                  CALL ctl_stop( "fld_init: yearly mean file must be in a yearly type of file: "//TRIM(sdjf%clrootname) )
388               ENDIF
389            ELSEIF( sdjf%nfreqh ==  -1 ) THEN   ! monthly mean
390               IF( sdjf%cltype == 'monthly' ) THEN            ! monthly file
391                  sdjf%nrec_a(1) = 1                                                       ! force to read the unique record
392                  llprevmth = .TRUE.                                                       ! use previous month file?
393                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
394               ELSE                                           ! yearly file
395                  sdjf%nrec_a(1) = 12                                                      ! force to read december mean
396                  llprevyr = .NOT. sdjf%ln_clim                                            ! use previous year  file?
397               ENDIF
398            ELSE                                ! higher frequency mean (in hours)
399               IF    ( sdjf%cltype      == 'monthly' ) THEN   ! monthly file
400                  sdjf%nrec_a(1) = NINT( 24 * nmonth_len(nmonth-1) / sdjf%nfreqh )         ! last record of previous month
401                  llprevmth = .TRUE.                                                       ! use previous month file?
402                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
403               ELSEIF( sdjf%cltype(1:4) == 'week'    ) THEN   ! weekly file
404                  llprevweek = .TRUE.                                                      ! use previous week  file?
405                  sdjf%nrec_a(1) = NINT( 24 * 7 / sdjf%nfreqh )                            ! last record of previous week
406                  isec_week = NINT(rday) * 7                                               ! add a shift toward previous week
407               ELSEIF( sdjf%cltype      == 'daily'   ) THEN   ! daily file
408                  sdjf%nrec_a(1) = NINT( 24 / sdjf%nfreqh )                                ! last record of previous day
409                  llprevday = .TRUE.                                                       ! use previous day   file?
410                  llprevmth = llprevday .AND. nday   == 1                                  ! use previous month file?
411                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
412               ELSE                                           ! yearly file
413                  sdjf%nrec_a(1) = NINT( 24 * nyear_len(0) / sdjf%nfreqh )                 ! last record of previous year
414                  llprevyr = .NOT. sdjf%ln_clim                                            ! use previous year  file?
415               ENDIF
416            ENDIF
417         ENDIF
418         !
419         IF ( sdjf%cltype(1:4) == 'week' ) THEN
420            isec_week = isec_week + ksec_week( sdjf%cltype(6:8) )   ! second since the beginning of the week
421            llprevmth = isec_week > nsec_month                      ! longer time since the beginning of the week than the month
422            llprevyr  = llprevmth .AND. nmonth == 1
423         ENDIF
424         llprev = llprevyr .OR. llprevmth .OR. llprevweek .OR. llprevday
425         !
426         iyear  = nyear  - COUNT((/llprevyr /))
427         imonth = nmonth - COUNT((/llprevmth/)) + 12 * COUNT((/llprevyr /))
428         iday   = nday   - COUNT((/llprevday/)) + nmonth_len(nmonth-1) * COUNT((/llprevmth/)) - isec_week / NINT(rday)
429         !
430         CALL fld_clopn( sdjf, iyear, imonth, iday, .NOT. llprev )
431         !
432         ! if previous year/month/day file does not exist, we switch to the current year/month/day
433         IF( llprev .AND. sdjf%num <= 0 ) THEN
434            CALL ctl_warn( 'previous year/month/week/day file: '//TRIM(sdjf%clrootname)//   &
435               &           ' not present -> back to current year/month/week/day' )
436            ! we force to read the first record of the current year/month/day instead of last record of previous year/month/day
437            llprev = .FALSE.
438            sdjf%nrec_a(1) = 1
439            CALL fld_clopn( sdjf )
440         ENDIF
441         !
442         IF( llprev ) THEN   ! check if the record sdjf%nrec_a(1) exists in the file
443            idvar = iom_varid( sdjf%num, sdjf%clvar )                                        ! id of the variable sdjf%clvar
444            IF( idvar <= 0 )   RETURN
445            inrec = iom_file( sdjf%num )%dimsz( iom_file( sdjf%num )%ndims(idvar), idvar )   ! size of the last dim of idvar
446            sdjf%nrec_a(1) = MIN( sdjf%nrec_a(1), inrec )   ! make sure we select an existing record
447         ENDIF
448         !
449         ! read before data in after arrays(as we will swap it later)
450         IF( PRESENT(jpk_bdy) ) THEN
451            CALL fld_get( sdjf, map, jpk_bdy, fvl )
452         ELSE
453            CALL fld_get( sdjf, map, ldxios = luxios)
454         ENDIF
455         !
456         clfmt = "('   fld_init : time-interpolation for ', a, ' read previous record = ', i6, ' at time = ', f7.2, ' days')"
457         IF(lwp) WRITE(numout, clfmt) TRIM(sdjf%clvar), sdjf%nrec_a(1), REAL(sdjf%nrec_a(2),wp)/rday
458         !
459      ENDIF
460      !
461   END SUBROUTINE fld_init
462
463
464   SUBROUTINE fld_rec( kn_fsbc, sdjf, ldbefore, kit, kt_offset )
465      !!---------------------------------------------------------------------
466      !!                    ***  ROUTINE fld_rec  ***
467      !!
468      !! ** Purpose : Compute
469      !!              if sdjf%ln_tint = .TRUE.
470      !!                  nrec_a: record number and its time (nrec_b is obtained from nrec_a when swapping)
471      !!              if sdjf%ln_tint = .FALSE.
472      !!                  nrec_a(1): record number
473      !!                  nrec_b(2) and nrec_a(2): time of the beginning and end of the record (for print only)
474      !!----------------------------------------------------------------------
475      INTEGER  , INTENT(in   )           ::   kn_fsbc   ! sbc computation period (in time step)
476      TYPE(FLD), INTENT(inout)           ::   sdjf      ! input field related variables
477      LOGICAL  , INTENT(in   ), OPTIONAL ::   ldbefore  ! sent back before record values (default = .FALSE.)
478      INTEGER  , INTENT(in   ), OPTIONAL ::   kit       ! index of barotropic subcycle
479      !                                                 ! used only if sdjf%ln_tint = .TRUE.
480      INTEGER  , INTENT(in   ), OPTIONAL ::   kt_offset ! Offset of required time level compared to "now"
481      !                                                 !   time level in units of time steps.
482      !
483      LOGICAL  ::   llbefore    ! local definition of ldbefore
484      INTEGER  ::   iendrec     ! end of this record (in seconds)
485      INTEGER  ::   imth        ! month number
486      INTEGER  ::   ifreq_sec   ! frequency mean (in seconds)
487      INTEGER  ::   isec_week   ! number of seconds since the start of the weekly file
488      INTEGER  ::   it_offset   ! local time offset variable
489      REAL(wp) ::   ztmp        ! temporary variable
490      !!----------------------------------------------------------------------
491      !
492      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
493      !
494      IF( PRESENT(ldbefore) ) THEN   ;   llbefore = ldbefore .AND. sdjf%ln_tint   ! needed only if sdjf%ln_tint = .TRUE.
495      ELSE                           ;   llbefore = .FALSE.
496      ENDIF
497      !
498      IF ( nn_components == jp_iam_sas ) THEN   ;   it_offset = nn_fsbc
499      ELSE                                      ;   it_offset = 0
500      ENDIF
501      IF( PRESENT(kt_offset) )   it_offset = kt_offset
502      IF( PRESENT(kit) ) THEN   ;   it_offset = ( kit + it_offset ) * NINT( rdt/REAL(nn_baro,wp) )
503      ELSE                      ;   it_offset =         it_offset   * NINT(       rdt            )
504      ENDIF
505      !
506      !                                      ! =========== !
507      IF    ( sdjf%nfreqh == -12 ) THEN      ! yearly mean
508         !                                   ! =========== !
509         !
510         IF( sdjf%ln_tint ) THEN                 ! time interpolation, shift by 1/2 record
511            !
512            !                  INT( ztmp )
513            !                     /|\
514            !                    1 |    *----
515            !                    0 |----(             
516            !                      |----+----|--> time
517            !                      0   /|\   1   (nday/nyear_len(1))
518            !                           |   
519            !                           |   
520            !       forcing record :    1
521            !                           
522            ztmp = REAL( nsec_year, wp ) / ( REAL( nyear_len(1), wp ) * rday ) + 0.5 &
523           &       + REAL( it_offset, wp ) / ( REAL( nyear_len(1), wp ) * rday )
524            sdjf%nrec_a(1) = 1 + INT( ztmp ) - COUNT((/llbefore/))
525            ! swap at the middle of the year
526            IF( llbefore ) THEN   ;   sdjf%nrec_a(2) = nsec1jan000 - (1 - INT(ztmp)) * NINT(0.5 * rday) * nyear_len(0) + &
527                                    & INT(ztmp) * NINT( 0.5 * rday) * nyear_len(1) 
528            ELSE                  ;   sdjf%nrec_a(2) = nsec1jan000 + (1 - INT(ztmp)) * NINT(0.5 * rday) * nyear_len(1) + &
529                                    & INT(ztmp) * INT(rday) * nyear_len(1) + INT(ztmp) * NINT( 0.5 * rday) * nyear_len(2) 
530            ENDIF
531         ELSE                                    ! no time interpolation
532            sdjf%nrec_a(1) = 1
533            sdjf%nrec_a(2) = NINT(rday) * nyear_len(1) + nsec1jan000   ! swap at the end    of the year
534            sdjf%nrec_b(2) = nsec1jan000                               ! beginning of the year (only for print)
535         ENDIF
536         !
537         !                                   ! ============ !
538      ELSEIF( sdjf%nfreqh ==  -1 ) THEN      ! monthly mean !
539         !                                   ! ============ !
540         !
541         IF( sdjf%ln_tint ) THEN                 ! time interpolation, shift by 1/2 record
542            !
543            !                  INT( ztmp )
544            !                     /|\
545            !                    1 |    *----
546            !                    0 |----(             
547            !                      |----+----|--> time
548            !                      0   /|\   1   (nday/nmonth_len(nmonth))
549            !                           |   
550            !                           |   
551            !       forcing record :  nmonth
552            !                           
553            ztmp = REAL( nsec_month, wp ) / ( REAL( nmonth_len(nmonth), wp ) * rday ) + 0.5 &
554           &       + REAL( it_offset, wp ) / ( REAL( nmonth_len(nmonth), wp ) * rday )
555            imth = nmonth + INT( ztmp ) - COUNT((/llbefore/))
556            IF( sdjf%cltype == 'monthly' ) THEN   ;   sdjf%nrec_a(1) = 1 + INT( ztmp ) - COUNT((/llbefore/))
557            ELSE                                  ;   sdjf%nrec_a(1) = imth
558            ENDIF
559            sdjf%nrec_a(2) = nmonth_half(   imth ) + nsec1jan000   ! swap at the middle of the month
560         ELSE                                    ! no time interpolation
561            IF( sdjf%cltype == 'monthly' ) THEN   ;   sdjf%nrec_a(1) = 1
562            ELSE                                  ;   sdjf%nrec_a(1) = nmonth
563            ENDIF
564            sdjf%nrec_a(2) =  nmonth_end(nmonth  ) + nsec1jan000   ! swap at the end    of the month
565            sdjf%nrec_b(2) =  nmonth_end(nmonth-1) + nsec1jan000   ! beginning of the month (only for print)
566         ENDIF
567         !
568         !                                   ! ================================ !
569      ELSE                                   ! higher frequency mean (in hours)
570         !                                   ! ================================ !
571         !
572         ifreq_sec = NINT( sdjf%nfreqh * 3600 )                                         ! frequency mean (in seconds)
573         IF( sdjf%cltype(1:4) == 'week' )   isec_week = ksec_week( sdjf%cltype(6:8) )   ! since the first day of the current week
574         ! number of second since the beginning of the file
575         IF(     sdjf%cltype      == 'monthly' ) THEN   ;   ztmp = REAL(nsec_month,wp)  ! since the first day of the current month
576         ELSEIF( sdjf%cltype(1:4) == 'week'    ) THEN   ;   ztmp = REAL(isec_week ,wp)  ! since the first day of the current week
577         ELSEIF( sdjf%cltype      == 'daily'   ) THEN   ;   ztmp = REAL(nsec_day  ,wp)  ! since 00h of the current day
578         ELSE                                           ;   ztmp = REAL(nsec_year ,wp)  ! since 00h on Jan 1 of the current year
579         ENDIF
580         ztmp = ztmp + 0.5 * REAL(kn_fsbc - 1, wp) * rdt + REAL( it_offset, wp )        ! centrered in the middle of sbc time step
581         ztmp = ztmp + 0.01 * rdt                                                       ! avoid truncation error
582         IF( sdjf%ln_tint ) THEN                ! time interpolation, shift by 1/2 record
583            !
584            !          INT( ztmp/ifreq_sec + 0.5 )
585            !                     /|\
586            !                    2 |        *-----(
587            !                    1 |  *-----(
588            !                    0 |--(             
589            !                      |--+--|--+--|--+--|--> time
590            !                      0 /|\ 1 /|\ 2 /|\ 3    (ztmp/ifreq_sec)
591            !                         |     |     |
592            !                         |     |     |
593            !       forcing record :  1     2     3
594            !                   
595            ztmp= ztmp / REAL(ifreq_sec, wp) + 0.5
596         ELSE                                   ! no time interpolation
597            !
598            !           INT( ztmp/ifreq_sec )
599            !                     /|\
600            !                    2 |           *-----(
601            !                    1 |     *-----(
602            !                    0 |-----(             
603            !                      |--+--|--+--|--+--|--> time
604            !                      0 /|\ 1 /|\ 2 /|\ 3    (ztmp/ifreq_sec)
605            !                         |     |     |
606            !                         |     |     |
607            !       forcing record :  1     2     3
608            !                           
609            ztmp= ztmp / REAL(ifreq_sec, wp)
610         ENDIF
611         sdjf%nrec_a(1) = 1 + INT( ztmp ) - COUNT((/llbefore/))   ! record number to be read
612
613         iendrec = ifreq_sec * sdjf%nrec_a(1) + nsec1jan000       ! end of this record (in second)
614         ! add the number of seconds between 00h Jan 1 and the end of previous month/week/day (ok if nmonth=1)
615         IF( sdjf%cltype      == 'monthly' )   iendrec = iendrec + NINT(rday) * SUM(nmonth_len(1:nmonth -1))
616         IF( sdjf%cltype(1:4) == 'week'    )   iendrec = iendrec + ( nsec_year - isec_week )
617         IF( sdjf%cltype      == 'daily'   )   iendrec = iendrec + NINT(rday) * ( nday_year - 1 )
618         IF( sdjf%ln_tint ) THEN
619             sdjf%nrec_a(2) = iendrec - ifreq_sec / 2        ! swap at the middle of the record
620         ELSE
621             sdjf%nrec_a(2) = iendrec                        ! swap at the end    of the record
622             sdjf%nrec_b(2) = iendrec - ifreq_sec            ! beginning of the record (only for print)
623         ENDIF
624         !
625      ENDIF
626      !
627   END SUBROUTINE fld_rec
628
629
630   SUBROUTINE fld_get( sdjf, map, jpk_bdy, fvl, ldxios)
631      !!---------------------------------------------------------------------
632      !!                    ***  ROUTINE fld_get  ***
633      !!
634      !! ** Purpose :   read the data
635      !!----------------------------------------------------------------------
636      TYPE(FLD)        , INTENT(inout) ::   sdjf   ! input field related variables
637      TYPE(MAP_POINTER), INTENT(in   ) ::   map    ! global-to-local mapping indices
638      INTEGER  , INTENT(in), OPTIONAL  ::   jpk_bdy ! number of vertical levels in the bdy data
639      LOGICAL  , INTENT(in), OPTIONAL  ::   fvl     ! number of vertical levels in the bdy data
640      LOGICAL  ,             OPTIONAL  ::   ldxios  ! xios I/O
641      !
642      INTEGER ::   ipk      ! number of vertical levels of sdjf%fdta ( 2D: ipk=1 ; 3D: ipk=jpk )
643      INTEGER ::   iw       ! index into wgts array
644      INTEGER ::   ipdom    ! index of the domain
645      INTEGER ::   idvar    ! variable ID
646      INTEGER ::   idmspc   ! number of spatial dimensions
647      LOGICAL ::   lmoor    ! C1D case: point data
648      LOGICAL ::   luxios   ! local variable for ldxios
649      REAL(wp)::   amaxval, aminval
650      !!---------------------------------------------------------------------
651      !
652
653      luxios = .FALSE.
654      IF( PRESENT(ldxios) ) luxios = ldxios
655      ipk = SIZE( sdjf%fnow, 3 )
656      !
657      IF( ASSOCIATED(map%ptr) ) THEN
658         IF( PRESENT(jpk_bdy) ) THEN
659            IF( sdjf%ln_tint ) THEN   ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fdta(:,:,:,2),                &
660                                                        sdjf%nrec_a(1), map, sdjf%igrd, sdjf%ibdy, jpk_bdy, fvl )
661            ELSE                      ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fnow(:,:,:  ),                &
662                                                        sdjf%nrec_a(1), map, sdjf%igrd, sdjf%ibdy, jpk_bdy, fvl )
663            ENDIF
664         ELSE
665            IF( sdjf%ln_tint ) THEN   ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fdta(:,:,:,2), sdjf%nrec_a(1), map )
666            ELSE                      ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fnow(:,:,:  ), sdjf%nrec_a(1), map )
667            ENDIF
668         ENDIF       
669      ELSE IF( LEN(TRIM(sdjf%wgtname)) > 0 ) THEN
670         CALL wgt_list( sdjf, iw, lxios_blkw )
671         IF( sdjf%ln_tint ) THEN   ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw, ipk, sdjf%fdta(:,:,:,2),          & 
672            &                                                                          sdjf%nrec_a(1), sdjf%lsmname )
673         ELSE                      ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw, ipk, sdjf%fnow(:,:,:  ),          &
674            &                                                                          sdjf%nrec_a(1), sdjf%lsmname )
675         ENDIF
676      ELSE
677         IF( SIZE(sdjf%fnow, 1) == jpi ) THEN   ;   ipdom = jpdom_data
678         ELSE                                   ;   ipdom = jpdom_unknown
679         ENDIF
680         ! C1D case: If product of spatial dimensions == ipk, then x,y are of
681         ! size 1 (point/mooring data): this must be read onto the central grid point
682         idvar  = iom_varid( sdjf%num, sdjf%clvar )
683         idmspc = iom_file ( sdjf%num )%ndims( idvar )
684         IF( iom_file( sdjf%num )%luld( idvar ) )   idmspc = idmspc - 1
685         lmoor  = (  idmspc == 0 .OR. PRODUCT( iom_file( sdjf%num )%dimsz( 1:MAX(idmspc,1) ,idvar ) ) == ipk  )
686         !
687         SELECT CASE( ipk )
688         CASE(1)
689            IF( lk_c1d .AND. lmoor ) THEN
690               IF( sdjf%ln_tint ) THEN
691                  CALL iom_get( sdjf%num, sdjf%clvar, sdjf%fdta(2,2,1,2), sdjf%nrec_a(1), ldxios = luxios )
692                  CALL lbc_lnk( 'fldread', sdjf%fdta(:,:,1,2),'Z',1. )
693               ELSE
694                  CALL iom_get( sdjf%num, sdjf%clvar, sdjf%fnow(2,2,1  ), sdjf%nrec_a(1), ldxios = luxios )
695                  CALL lbc_lnk( 'fldread', sdjf%fnow(:,:,1  ),'Z',1. )
696               ENDIF
697            ELSE
698               IF( sdjf%ln_tint ) THEN   ;   CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fdta(:,:,1,2), sdjf%nrec_a(1), ldxios = luxios )
699               ELSE                      ;   CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fnow(:,:,1  ), sdjf%nrec_a(1), ldxios = luxios )
700               ENDIF
701            ENDIF
702         CASE DEFAULT
703            IF (lk_c1d .AND. lmoor ) THEN
704               IF( sdjf%ln_tint ) THEN
705                  CALL iom_get( sdjf%num, jpdom_unknown, sdjf%clvar, sdjf%fdta(2,2,:,2), sdjf%nrec_a(1), ldxios = luxios)
706                  CALL lbc_lnk( 'fldread', sdjf%fdta(:,:,:,2),'Z',1. )
707               ELSE
708                  CALL iom_get( sdjf%num, jpdom_unknown, sdjf%clvar, sdjf%fnow(2,2,:  ), sdjf%nrec_a(1), ldxios = luxios)
709                  CALL lbc_lnk( 'fldread', sdjf%fnow(:,:,:  ),'Z',1. )
710               ENDIF
711            ELSE
712               IF( sdjf%ln_tint ) THEN   ;   CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fdta(:,:,:,2), sdjf%nrec_a(1), ldxios = luxios )
713               ELSE                      ;   CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fnow(:,:,:  ), sdjf%nrec_a(1), ldxios = luxios )
714               ENDIF
715            ENDIF
716         END SELECT
717      ENDIF
718      !
719      sdjf%rotn(2) = .false.   ! vector not yet rotated
720      !
721   END SUBROUTINE fld_get
722
723   SUBROUTINE fld_map( num, clvar, dta, nrec, map, igrd, ibdy, jpk_bdy, fvl )
724      !!---------------------------------------------------------------------
725      !!                    ***  ROUTINE fld_map  ***
726      !!
727      !! ** Purpose :   read global data from file and map onto local data
728      !!                using a general mapping (for open boundaries)
729      !!----------------------------------------------------------------------
730
731      USE bdy_oce, ONLY: ln_bdy, idx_bdy, dta_global, dta_global_z, dta_global_dz, dta_global2, dta_global2_z, dta_global2_dz                 ! workspace to read in global data arrays
732
733      INTEGER                   , INTENT(in ) ::   num     ! stream number
734      CHARACTER(LEN=*)          , INTENT(in ) ::   clvar   ! variable name
735      REAL(wp), DIMENSION(:,:,:), INTENT(out) ::   dta     ! output field on model grid (2 dimensional)
736      INTEGER                   , INTENT(in ) ::   nrec    ! record number to read (ie time slice)
737      TYPE(MAP_POINTER)         , INTENT(in ) ::   map     ! global-to-local mapping indices
738      INTEGER  , INTENT(in), OPTIONAL         ::   igrd, ibdy, jpk_bdy  ! grid type, set number and number of vertical levels in the bdy data
739      LOGICAL  , INTENT(in), OPTIONAL         ::   fvl     ! grid type, set number and number of vertical levels in the bdy data
740      INTEGER                                 ::   jpkm1_bdy! number of vertical levels in the bdy data minus 1
741      !!
742      INTEGER                                 ::   ipi      ! length of boundary data on local process
743      INTEGER                                 ::   ipj      ! length of dummy dimension ( = 1 )
744      INTEGER                                 ::   ipk      ! number of vertical levels of dta ( 2D: ipk=1 ; 3D: ipk=jpk )
745      INTEGER                                 ::   ilendta  ! length of data in file
746      INTEGER                                 ::   idvar    ! variable ID
747      INTEGER                                 ::   ib, ik, ji, jj   ! loop counters
748      INTEGER                                 ::   ierr
749      REAL(wp)                                ::   fv          ! fillvalue
750      REAL(wp), POINTER, DIMENSION(:,:,:)     ::   dta_read    ! work space for global data
751      REAL(wp), POINTER, DIMENSION(:,:,:)     ::   dta_read_z  ! work space for global data
752      REAL(wp), POINTER, DIMENSION(:,:,:)     ::   dta_read_dz ! work space for global data
753      !!---------------------------------------------------------------------
754      !
755      ipi = SIZE( dta, 1 )
756      ipj = 1
757      ipk = SIZE( dta, 3 )
758      !
759      idvar   = iom_varid( num, clvar )
760      ilendta = iom_file(num)%dimsz(1,idvar)
761
762      IF ( ln_bdy ) THEN
763         ipj = iom_file(num)%dimsz(2,idvar)
764         IF( map%ll_unstruc) THEN   ! unstructured open boundary data file
765            dta_read => dta_global
766            IF( PRESENT(jpk_bdy) ) THEN
767               IF( jpk_bdy>0 ) THEN
768                  dta_read_z => dta_global_z
769                  dta_read_dz => dta_global_dz
770                  jpkm1_bdy = jpk_bdy-1
771               ENDIF
772            ENDIF
773         ELSE                       ! structured open boundary file
774            dta_read => dta_global2
775            IF( PRESENT(jpk_bdy) ) THEN
776               IF( jpk_bdy>0 ) THEN
777                  dta_read_z => dta_global2_z
778                  dta_read_dz => dta_global2_dz
779                  jpkm1_bdy = jpk_bdy-1
780               ENDIF
781            ENDIF
782         ENDIF
783      ENDIF
784
785      IF(lwp) WRITE(numout,*) 'Dim size for ',        TRIM(clvar),' is ', ilendta
786      IF(lwp) WRITE(numout,*) 'Number of levels for ',TRIM(clvar),' is ', ipk
787      !
788      SELECT CASE( ipk )
789      CASE(1)        ;   
790      CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1    ), nrec )
791         IF ( map%ll_unstruc) THEN ! unstructured open boundary data file
792            DO ib = 1, ipi
793              DO ik = 1, ipk
794                dta(ib,1,ik) =  dta_read(map%ptr(ib),1,ik)
795              END DO
796            END DO
797         ELSE ! we assume that this is a structured open boundary file
798            DO ib = 1, ipi
799               jj=1+floor(REAL(map%ptr(ib)-1)/REAL(ilendta))
800               ji=map%ptr(ib)-(jj-1)*ilendta
801               DO ik = 1, ipk
802                  dta(ib,1,ik) =  dta_read(ji,jj,ik)
803               END DO
804            END DO
805         ENDIF
806
807      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
808      ! Do we include something here to adjust barotropic velocities !
809      ! in case of a depth difference between bdy files and          !
810      ! bathymetry in the case ln_full_vel = .false. and jpk_bdy>0?  !
811      ! [as the enveloping and parital cells could change H]         !
812      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
813
814      CASE DEFAULT   ;
815
816      IF( PRESENT(jpk_bdy) .AND. jpk_bdy>0 ) THEN       ! boundary data not on model grid: vertical interpolation
817         CALL iom_getatt(num, '_FillValue', fv, cdvar=clvar )
818         dta_read(:,:,:) = -ABS(fv)
819         dta_read_z(:,:,:) = 0._wp
820         dta_read_dz(:,:,:) = 0._wp
821         CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1:jpk_bdy), nrec )
822         SELECT CASE( igrd )                 
823            CASE(1)
824               CALL iom_get ( num, jpdom_unknown, 'gdept', dta_read_z(1:ilendta,1:ipj,1:jpk_bdy) )
825               CALL iom_get ( num, jpdom_unknown, 'e3t',  dta_read_dz(1:ilendta,1:ipj,1:jpk_bdy) )
826            CASE(2) 
827               CALL iom_get ( num, jpdom_unknown, 'gdepu', dta_read_z(1:ilendta,1:ipj,1:jpk_bdy) )
828               CALL iom_get ( num, jpdom_unknown, 'e3u',  dta_read_dz(1:ilendta,1:ipj,1:jpk_bdy) )
829            CASE(3)
830               CALL iom_get ( num, jpdom_unknown, 'gdepv', dta_read_z(1:ilendta,1:ipj,1:jpk_bdy) )
831               CALL iom_get ( num, jpdom_unknown, 'e3v',  dta_read_dz(1:ilendta,1:ipj,1:jpk_bdy) )
832         END SELECT
833
834      IF ( ln_bdy ) & 
835         CALL fld_bdy_interp(dta_read, dta_read_z, dta_read_dz, map, jpk_bdy, igrd, ibdy, fv, dta, fvl, ilendta)
836
837      ELSE ! boundary data assumed to be on model grid
838         
839         CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1:ipk), nrec )                   
840         IF ( map%ll_unstruc) THEN ! unstructured open boundary data file
841            DO ib = 1, ipi
842              DO ik = 1, ipk
843                dta(ib,1,ik) =  dta_read(map%ptr(ib),1,ik)
844              END DO
845            END DO
846         ELSE ! we assume that this is a structured open boundary file
847            DO ib = 1, ipi
848               jj=1+floor(REAL(map%ptr(ib)-1)/REAL(ilendta))
849               ji=map%ptr(ib)-(jj-1)*ilendta
850               DO ik = 1, ipk
851                  dta(ib,1,ik) =  dta_read(ji,jj,ik)
852               END DO
853            END DO
854         ENDIF
855      ENDIF ! PRESENT(jpk_bdy)
856      END SELECT
857
858   END SUBROUTINE fld_map
859   
860   SUBROUTINE fld_bdy_interp(dta_read, dta_read_z, dta_read_dz, map, jpk_bdy, igrd, ibdy, fv, dta, fvl, ilendta)
861
862      !!---------------------------------------------------------------------
863      !!                    ***  ROUTINE fld_bdy_interp  ***
864      !!
865      !! ** Purpose :   on the fly vertical interpolation to allow the use of
866      !!                boundary data from non-native vertical grid
867      !!----------------------------------------------------------------------
868      USE bdy_oce, ONLY:  idx_bdy         ! indexing for map <-> ij transformation
869
870      REAL(wp), POINTER, DIMENSION(:,:,:), INTENT(in )     ::   dta_read      ! work space for global data
871      REAL(wp), POINTER, DIMENSION(:,:,:), INTENT(in )     ::   dta_read_z    ! work space for global data
872      REAL(wp), POINTER, DIMENSION(:,:,:), INTENT(in )     ::   dta_read_dz   ! work space for global data
873      REAL(wp) , INTENT(in)                                ::   fv            ! fillvalue and alternative -ABS(fv)
874      REAL(wp), DIMENSION(:,:,:), INTENT(out) ::   dta                        ! output field on model grid (2 dimensional)
875      TYPE(MAP_POINTER)         , INTENT(in ) ::   map                        ! global-to-local mapping indices
876      LOGICAL  , INTENT(in), OPTIONAL         ::   fvl                        ! grid type, set number and number of vertical levels in the bdy data
877      INTEGER  , INTENT(in)                   ::   igrd, ibdy, jpk_bdy        ! number of levels in bdy data
878      INTEGER  , INTENT(in)                   ::   ilendta                    ! length of data in file
879      !!
880      INTEGER                                 ::   ipi                        ! length of boundary data on local process
881      INTEGER                                 ::   ipj                        ! length of dummy dimension ( = 1 )
882      INTEGER                                 ::   ipk                        ! number of vertical levels of dta ( 2D: ipk=1 ; 3D: ipk=jpk )
883      INTEGER                                 ::   jpkm1_bdy                  ! number of levels in bdy data minus 1
884      INTEGER                                 ::   ib, ik, ikk                ! loop counters
885      INTEGER                                 ::   ji, jj, zij, zjj           ! temporary indices
886      REAL(wp)                                ::   zl, zi, zh                 ! tmp variable for current depth and interpolation factor
887      REAL(wp)                                ::   fv_alt, ztrans, ztrans_new ! fillvalue and alternative -ABS(fv)
888      CHARACTER (LEN=10)                      ::   ibstr
889      !!---------------------------------------------------------------------
890     
891
892      ipi       = SIZE( dta, 1 )
893      ipj       = SIZE( dta_read, 2 )
894      ipk       = SIZE( dta, 3 )
895      jpkm1_bdy = jpk_bdy-1
896     
897      fv_alt = -ABS(fv)  ! set _FillValue < 0 as we make use of MAXVAL and MAXLOC later
898      DO ib = 1, ipi
899            zij = idx_bdy(ibdy)%nbi(ib,igrd)
900            zjj = idx_bdy(ibdy)%nbj(ib,igrd)
901         IF(narea==2) WRITE(*,*) 'MAPI', ib, igrd, map%ptr(ib), narea-1, zij, zjj
902      ENDDO
903      !
904      IF ( map%ll_unstruc ) THEN                            ! unstructured open boundary data file
905
906         DO ib = 1, ipi
907            DO ik = 1, jpk_bdy
908               IF( ( dta_read(map%ptr(ib),1,ik) == fv ) ) THEN
909                  dta_read_z(map%ptr(ib),1,ik)  = fv_alt ! safety: put fillvalue into external depth field so consistent with data
910                  dta_read_dz(map%ptr(ib),1,ik) = 0._wp  ! safety: put 0._wp into external thickness factors to ensure transport is correct
911               ENDIF
912            ENDDO
913         ENDDO 
914
915         DO ib = 1, ipi
916            zij = idx_bdy(ibdy)%nbi(ib,igrd)
917            zjj = idx_bdy(ibdy)%nbj(ib,igrd)
918            zh  = SUM(dta_read_dz(map%ptr(ib),1,:) )
919            ! Warnings to flag differences in the input and model topgraphy - is this useful/necessary?
920            SELECT CASE( igrd )                         
921               CASE(1)
922                  IF( ABS( (zh - ht_n(zij,zjj)) / ht_n(zij,zjj)) * tmask(zij,zjj,1) > 0.01_wp ) THEN
923                     WRITE(ibstr,"(I10.10)") map%ptr(ib) 
924                     CALL ctl_warn('fld_bdy_interp: T depths differ between grids at BDY point '//TRIM(ibstr)//' by more than 1%')
925                 !   IF(lwp) WRITE(*,*) 'DEPTHT', zh, sum(e3t_n(zij,zjj,:), mask=tmask(zij,zjj,:)==1),  ht_n(zij,zjj), map%ptr(ib), ib, zij, zjj
926                  ENDIF
927               CASE(2)
928                  IF( ABS( (zh - hu_n(zij,zjj)) * r1_hu_n(zij,zjj)) * umask(zij,zjj,1) > 0.01_wp ) THEN
929                     WRITE(ibstr,"(I10.10)") map%ptr(ib) 
930                     CALL ctl_warn('fld_bdy_interp: U depths differ between grids at BDY point '//TRIM(ibstr)//' by more than 1%')
931                     IF(lwp) WRITE(*,*) 'DEPTHU', zh, sum(e3u_n(zij,zjj,:), mask=umask(zij,zjj,:)==1),  sum(umask(zij,zjj,:)), &
932                       &                hu_n(zij,zjj), map%ptr(ib), ib, zij, zjj, narea-1  , &
933                        &                dta_read(map%ptr(ib),1,:)
934                  ENDIF
935               CASE(3)
936                  IF( ABS( (zh - hv_n(zij,zjj)) * r1_hv_n(zij,zjj)) * vmask(zij,zjj,1) > 0.01_wp ) THEN
937                     WRITE(ibstr,"(I10.10)") map%ptr(ib) 
938                     CALL ctl_warn('fld_bdy_interp: V depths differ between grids at BDY point '//TRIM(ibstr)//' by more than 1%')
939                  ENDIF
940            END SELECT
941            DO ik = 1, ipk                     
942               SELECT CASE( igrd )                       
943                  CASE(1)
944                     zl =  gdept_n(zij,zjj,ik)                                          ! if using in step could use fsdept instead of gdept_n?
945                  CASE(2)
946                     IF(ln_sco) THEN
947                        zl =  ( gdept_n(zij,zjj,ik) + gdept_n(zij+1,zjj,ik) ) * 0.5_wp  ! if using in step could use fsdept instead of gdept_n?
948                     ELSE
949                        zl =  MIN( gdept_n(zij,zjj,ik), gdept_n(zij+1,zjj,ik) ) 
950                     ENDIF
951                  CASE(3)
952                     IF(ln_sco) THEN
953                        zl =  ( gdept_n(zij,zjj,ik) + gdept_n(zij,zjj+1,ik) ) * 0.5_wp  ! if using in step could use fsdept instead of gdept_n?
954                     ELSE
955                        zl =  MIN( gdept_n(zij,zjj,ik), gdept_n(zij,zjj+1,ik) )
956                     ENDIF
957               END SELECT
958               IF( zl < dta_read_z(map%ptr(ib),1,1) ) THEN                                         ! above the first level of external data
959                  dta(ib,1,ik) =  dta_read(map%ptr(ib),1,1)
960               ELSEIF( zl > MAXVAL(dta_read_z(map%ptr(ib),1,:),1) ) THEN                           ! below the last level of external data
961                  dta(ib,1,ik) =  dta_read(map%ptr(ib),1,MAXLOC(dta_read_z(map%ptr(ib),1,:),1))
962               ELSE                                                                                ! inbetween : vertical interpolation between ikk & ikk+1
963                  DO ikk = 1, jpkm1_bdy                                                            ! when  gdept_n(ikk) < zl < gdept_n(ikk+1)
964                     IF( ( (zl-dta_read_z(map%ptr(ib),1,ikk)) * (zl-dta_read_z(map%ptr(ib),1,ikk+1)) <= 0._wp) &
965                    &    .AND. (dta_read_z(map%ptr(ib),1,ikk+1) /= fv_alt)) THEN
966                        zi           = ( zl - dta_read_z(map%ptr(ib),1,ikk) ) / &
967                       &               ( dta_read_z(map%ptr(ib),1,ikk+1) - dta_read_z(map%ptr(ib),1,ikk) )
968                        dta(ib,1,ik) = dta_read(map%ptr(ib),1,ikk) + &
969                       &               ( dta_read(map%ptr(ib),1,ikk+1) - dta_read(map%ptr(ib),1,ikk) ) * zi
970                     ENDIF
971                  END DO
972               ENDIF
973            END DO
974         END DO
975
976         IF(igrd == 2) THEN                                 ! do we need to adjust the transport term?
977            DO ib = 1, ipi
978              zij = idx_bdy(ibdy)%nbi(ib,igrd)
979              zjj = idx_bdy(ibdy)%nbj(ib,igrd)
980              zh  = SUM(dta_read_dz(map%ptr(ib),1,:) )
981              ztrans = 0._wp
982              ztrans_new = 0._wp
983              DO ik = 1, jpk_bdy                            ! calculate transport on input grid
984                  ztrans     = ztrans     + dta_read(map%ptr(ib),1,ik) * dta_read_dz(map%ptr(ib),1,ik)
985              ENDDO
986              DO ik = 1, ipk                                ! calculate transport on model grid
987                  ztrans_new = ztrans_new + dta(ib,1,ik) * e3u_n(zij,zjj,ik) * umask(zij,zjj,ik)
988              ENDDO
989              DO ik = 1, ipk                                ! make transport correction
990                 IF(fvl) THEN ! bdy data are total velocity so adjust bt transport term to match input data
991                    dta(ib,1,ik) = ( dta(ib,1,ik) + ( ztrans - ztrans_new ) * r1_hu_n(zij,zjj) ) * umask(zij,zjj,ik)
992                 ELSE ! we're just dealing with bc velocity so bt transport term should sum to zero
993                    IF( ABS(ztrans * r1_hu_n(zij,zjj)) > 0.01_wp ) &
994                   &   CALL ctl_warn('fld_bdy_interp: barotropic component of > 0.01 ms-1 found in baroclinic velocities at')
995                    dta(ib,1,ik) = dta(ib,1,ik) + ( 0._wp - ztrans_new ) * r1_hu_n(zij,zjj) * umask(zij,zjj,ik)
996                 ENDIF
997              ENDDO
998            ENDDO
999         ENDIF
1000
1001         IF(igrd == 3) THEN                                 ! do we need to adjust the transport term?
1002            DO ib = 1, ipi
1003              zij = idx_bdy(ibdy)%nbi(ib,igrd)
1004              zjj = idx_bdy(ibdy)%nbj(ib,igrd)
1005              zh  = SUM(dta_read_dz(map%ptr(ib),1,:) )
1006              ztrans = 0._wp
1007              ztrans_new = 0._wp
1008              DO ik = 1, jpk_bdy                            ! calculate transport on input grid
1009                  ztrans     = ztrans     + dta_read(map%ptr(ib),1,ik) * dta_read_dz(map%ptr(ib),1,ik)
1010              ENDDO
1011              DO ik = 1, ipk                                ! calculate transport on model grid
1012                  ztrans_new = ztrans_new + dta(ib,1,ik) * e3v_n(zij,zjj,ik) * vmask(zij,zjj,ik)
1013              ENDDO
1014              DO ik = 1, ipk                                ! make transport correction
1015                 IF(fvl) THEN ! bdy data are total velocity so adjust bt transport term to match input data
1016                    dta(ib,1,ik) = ( dta(ib,1,ik) + ( ztrans - ztrans_new ) * r1_hv_n(zij,zjj) ) * vmask(zij,zjj,ik)
1017                 ELSE ! we're just dealing with bc velocity so bt transport term should sum to zero
1018                    dta(ib,1,ik) = dta(ib,1,ik) + ( 0._wp - ztrans_new ) * r1_hv_n(zij,zjj) * vmask(zij,zjj,ik)
1019                 ENDIF
1020              ENDDO
1021            ENDDO
1022         ENDIF
1023 
1024      ELSE ! structured open boundary file
1025
1026         DO ib = 1, ipi
1027            jj=1+floor(REAL(map%ptr(ib)-1)/REAL(ilendta))
1028            ji=map%ptr(ib)-(jj-1)*ilendta
1029            DO ik = 1, jpk_bdy                     
1030               IF( ( dta_read(ji,jj,ik) == fv ) ) THEN
1031                  dta_read_z(ji,jj,ik)  = fv_alt ! safety: put fillvalue into external depth field so consistent with data
1032                  dta_read_dz(ji,jj,ik) = 0._wp  ! safety: put 0._wp into external thickness factors to ensure transport is correct
1033               ENDIF
1034            ENDDO
1035         ENDDO 
1036       
1037
1038         DO ib = 1, ipi
1039            jj=1+floor(REAL(map%ptr(ib)-1)/REAL(ilendta))
1040            ji=map%ptr(ib)-(jj-1)*ilendta
1041            zij = idx_bdy(ibdy)%nbi(ib,igrd)
1042            zjj = idx_bdy(ibdy)%nbj(ib,igrd)
1043            zh  = SUM(dta_read_dz(ji,jj,:) )
1044            ! Warnings to flag differences in the input and model topgraphy - is this useful/necessary?
1045            SELECT CASE( igrd )                         
1046               CASE(1)
1047                  IF( ABS( (zh - ht_n(zij,zjj)) / ht_n(zij,zjj)) * tmask(zij,zjj,1) > 0.01_wp ) THEN
1048                     WRITE(ibstr,"(I10.10)") map%ptr(ib) 
1049                     CALL ctl_warn('fld_bdy_interp: T depths differ between grids at BDY point '//TRIM(ibstr)//' by more than 1%')
1050                 !   IF(lwp) WRITE(*,*) 'DEPTHT', zh, sum(e3t_n(zij,zjj,:), mask=tmask(zij,zjj,:)==1),  ht_n(zij,zjj), map%ptr(ib), ib, zij, zjj
1051                  ENDIF
1052               CASE(2)
1053                  IF( ABS( (zh - hu_n(zij,zjj)) * r1_hu_n(zij,zjj)) * umask(zij,zjj,1) > 0.01_wp ) THEN
1054                     WRITE(ibstr,"(I10.10)") map%ptr(ib) 
1055                     CALL ctl_warn('fld_bdy_interp: U depths differ between grids at BDY point '//TRIM(ibstr)//' by more than 1%')
1056                  ENDIF
1057               CASE(3)
1058                  IF( ABS( (zh - hv_n(zij,zjj)) * r1_hv_n(zij,zjj)) * vmask(zij,zjj,1) > 0.01_wp ) THEN
1059                     WRITE(ibstr,"(I10.10)") map%ptr(ib) 
1060                     CALL ctl_warn('fld_bdy_interp: V depths differ between grids at BDY point '//TRIM(ibstr)//' by more than 1%')
1061                  ENDIF
1062            END SELECT
1063            DO ik = 1, ipk                     
1064               SELECT CASE( igrd )                          ! coded for sco - need zco and zps option using min
1065                  CASE(1)
1066                     zl =  gdept_n(zij,zjj,ik)              ! if using in step could use fsdept instead of gdept_n?
1067                  CASE(2)
1068                     IF(ln_sco) THEN
1069                        zl =  ( gdept_n(zij,zjj,ik) + gdept_n(zij+1,zjj,ik) ) * 0.5_wp  ! if using in step could use fsdept instead of gdept_n?
1070                     ELSE
1071                        zl =  MIN( gdept_n(zij,zjj,ik), gdept_n(zij+1,zjj,ik) )
1072                     ENDIF
1073                  CASE(3)
1074                     IF(ln_sco) THEN
1075                        zl =  ( gdept_n(zij,zjj,ik) + gdept_n(zij,zjj+1,ik) ) * 0.5_wp  ! if using in step could use fsdept instead of gdept_n?
1076                     ELSE
1077                        zl =  MIN( gdept_n(zij,zjj,ik), gdept_n(zij,zjj+1,ik) )
1078                     ENDIF
1079               END SELECT
1080               IF( zl < dta_read_z(ji,jj,1) ) THEN                                      ! above the first level of external data
1081                  dta(ib,1,ik) =  dta_read(ji,jj,1)
1082               ELSEIF( zl > MAXVAL(dta_read_z(ji,jj,:),1) ) THEN                        ! below the last level of external data
1083                  dta(ib,1,ik) =  dta_read(ji,jj,MAXLOC(dta_read_z(ji,jj,:),1))
1084               ELSE                                                                     ! inbetween : vertical interpolation between ikk & ikk+1
1085                  DO ikk = 1, jpkm1_bdy                                                 ! when  gdept_n(ikk) < zl < gdept_n(ikk+1)
1086                     IF( ( (zl-dta_read_z(ji,jj,ikk)) * (zl-dta_read_z(ji,jj,ikk+1)) <= 0._wp) &
1087                    &    .AND. (dta_read_z(ji,jj,ikk+1) /= fv_alt)) THEN
1088                        zi           = ( zl - dta_read_z(ji,jj,ikk) ) / &
1089                       &               ( dta_read_z(ji,jj,ikk+1) - dta_read_z(ji,jj,ikk) )
1090                        dta(ib,1,ik) = dta_read(ji,jj,ikk) + &
1091                       &               ( dta_read(ji,jj,ikk+1) - dta_read(ji,jj,ikk) ) * zi
1092                     ENDIF
1093                  END DO
1094               ENDIF
1095            END DO
1096         END DO
1097
1098         IF(igrd == 2) THEN                                 ! do we need to adjust the transport term?
1099            DO ib = 1, ipi
1100               jj=1+floor(REAL(map%ptr(ib)-1)/REAL(ilendta))
1101               ji=map%ptr(ib)-(jj-1)*ilendta
1102               zij = idx_bdy(ibdy)%nbi(ib,igrd)
1103               zjj = idx_bdy(ibdy)%nbj(ib,igrd)
1104               zh = SUM(dta_read_dz(ji,jj,:) )
1105               ztrans = 0._wp
1106               ztrans_new = 0._wp
1107               DO ik = 1, jpk_bdy                            ! calculate transport on input grid
1108                  ztrans = ztrans + dta_read(ji,jj,ik) * dta_read_dz(ji,jj,ik)
1109               ENDDO
1110               DO ik = 1, ipk                                ! calculate transport on model grid
1111                  ztrans_new = ztrans_new + dta(ib,1,ik) * e3u_n(zij,zjj,ik) * umask(zij,zjj,ik)
1112               ENDDO
1113               DO ik = 1, ipk                                ! make transport correction
1114                  IF(fvl) THEN ! bdy data are total velocity so adjust bt transport term to match input data
1115                     dta(ib,1,ik) = ( dta(ib,1,ik) + ( ztrans - ztrans_new ) * r1_hu_n(zij,zjj) ) * umask(zij,zjj,ik)
1116                  ELSE ! we're just dealing with bc velocity so bt transport term should sum to zero
1117                     dta(ib,1,ik) = ( dta(ib,1,ik) + ( 0._wp  - ztrans_new ) * r1_hu_n(zij,zjj) ) * umask(zij,zjj,ik)
1118                  ENDIF
1119               ENDDO
1120            ENDDO
1121         ENDIF
1122
1123         IF(igrd == 3) THEN                                 ! do we need to adjust the transport term?
1124            DO ib = 1, ipi
1125               jj  = 1+floor(REAL(map%ptr(ib)-1)/REAL(ilendta))
1126               ji  = map%ptr(ib)-(jj-1)*ilendta
1127               zij = idx_bdy(ibdy)%nbi(ib,igrd)
1128               zjj = idx_bdy(ibdy)%nbj(ib,igrd)
1129               zh  = SUM(dta_read_dz(ji,jj,:) )
1130               ztrans = 0._wp
1131               ztrans_new = 0._wp
1132               DO ik = 1, jpk_bdy                            ! calculate transport on input grid
1133                  ztrans     = ztrans     + dta_read(ji,jj,ik) * dta_read_dz(ji,jj,ik)
1134               ENDDO
1135               DO ik = 1, ipk                                ! calculate transport on model grid
1136                  ztrans_new = ztrans_new + dta(ib,1,ik) * e3v_n(zij,zjj,ik) * vmask(zij,zjj,ik)
1137               ENDDO
1138               DO ik = 1, ipk                                ! make transport correction
1139                  IF(fvl) THEN ! bdy data are total velocity so adjust bt transport term to match input data
1140                     dta(ib,1,ik) = ( dta(ib,1,ik) + ( ztrans - ztrans_new ) * r1_hv_n(zij,zjj) ) * vmask(zij,zjj,ik)
1141                  ELSE ! we're just dealing with bc velocity so bt transport term should sum to zero
1142                     dta(ib,1,ik) = ( dta(ib,1,ik) + ( 0._wp  - ztrans_new ) * r1_hv_n(zij,zjj) ) * vmask(zij,zjj,ik)
1143                  ENDIF
1144               ENDDO
1145            ENDDO
1146         ENDIF
1147
1148      ENDIF ! endif unstructured or structured
1149
1150   END SUBROUTINE fld_bdy_interp
1151
1152
1153   SUBROUTINE fld_rot( kt, sd )
1154      !!---------------------------------------------------------------------
1155      !!                    ***  ROUTINE fld_rot  ***
1156      !!
1157      !! ** Purpose :   Vector fields may need to be rotated onto the local grid direction
1158      !!----------------------------------------------------------------------
1159      INTEGER                , INTENT(in   ) ::   kt   ! ocean time step
1160      TYPE(FLD), DIMENSION(:), INTENT(inout) ::   sd   ! input field related variables
1161      !
1162      INTEGER ::   ju, jv, jk, jn  ! loop indices
1163      INTEGER ::   imf             ! size of the structure sd
1164      INTEGER ::   ill             ! character length
1165      INTEGER ::   iv              ! indice of V component
1166      CHARACTER (LEN=100)          ::   clcomp       ! dummy weight name
1167      REAL(wp), DIMENSION(jpi,jpj) ::   utmp, vtmp   ! temporary arrays for vector rotation
1168      !!---------------------------------------------------------------------
1169      !
1170      !! (sga: following code should be modified so that pairs arent searched for each time
1171      !
1172      imf = SIZE( sd )
1173      DO ju = 1, imf
1174         ill = LEN_TRIM( sd(ju)%vcomp )
1175         DO jn = 2-COUNT((/sd(ju)%ln_tint/)), 2
1176            IF( ill > 0 .AND. .NOT. sd(ju)%rotn(jn) ) THEN   ! find vector rotations required             
1177               IF( sd(ju)%vcomp(1:1) == 'U' ) THEN      ! east-west component has symbolic name starting with 'U'
1178                  ! look for the north-south component which has same symbolic name but with 'U' replaced with 'V'
1179                  clcomp = 'V' // sd(ju)%vcomp(2:ill)   ! works even if ill == 1
1180                  iv = -1
1181                  DO jv = 1, imf
1182                     IF( TRIM(sd(jv)%vcomp) == TRIM(clcomp) )   iv = jv
1183                  END DO
1184                  IF( iv > 0 ) THEN   ! fields ju and iv are two components which need to be rotated together
1185                     DO jk = 1, SIZE( sd(ju)%fnow, 3 )
1186                        IF( sd(ju)%ln_tint )THEN
1187                           CALL rot_rep( sd(ju)%fdta(:,:,jk,jn), sd(iv)%fdta(:,:,jk,jn), 'T', 'en->i', utmp(:,:) )
1188                           CALL rot_rep( sd(ju)%fdta(:,:,jk,jn), sd(iv)%fdta(:,:,jk,jn), 'T', 'en->j', vtmp(:,:) )
1189                           sd(ju)%fdta(:,:,jk,jn) = utmp(:,:)   ;   sd(iv)%fdta(:,:,jk,jn) = vtmp(:,:)
1190                        ELSE
1191                           CALL rot_rep( sd(ju)%fnow(:,:,jk  ), sd(iv)%fnow(:,:,jk  ), 'T', 'en->i', utmp(:,:) )
1192                           CALL rot_rep( sd(ju)%fnow(:,:,jk  ), sd(iv)%fnow(:,:,jk  ), 'T', 'en->j', vtmp(:,:) )
1193                           sd(ju)%fnow(:,:,jk   ) = utmp(:,:)   ;   sd(iv)%fnow(:,:,jk   ) = vtmp(:,:)
1194                        ENDIF
1195                     END DO
1196                     sd(ju)%rotn(jn) = .TRUE.               ! vector was rotated
1197                     IF( lwp .AND. kt == nit000 )   WRITE(numout,*)   &
1198                        &   'fld_read: vector pair ('//TRIM(sd(ju)%clvar)//', '//TRIM(sd(iv)%clvar)//') rotated on to model grid'
1199                  ENDIF
1200               ENDIF
1201            ENDIF
1202         END DO
1203       END DO
1204      !
1205   END SUBROUTINE fld_rot
1206
1207
1208   SUBROUTINE fld_clopn( sdjf, kyear, kmonth, kday, ldstop )
1209      !!---------------------------------------------------------------------
1210      !!                    ***  ROUTINE fld_clopn  ***
1211      !!
1212      !! ** Purpose :   update the file name and close/open the files
1213      !!----------------------------------------------------------------------
1214      TYPE(FLD)        , INTENT(inout) ::   sdjf     ! input field related variables
1215      INTEGER, OPTIONAL, INTENT(in   ) ::   kyear    ! year value
1216      INTEGER, OPTIONAL, INTENT(in   ) ::   kmonth   ! month value
1217      INTEGER, OPTIONAL, INTENT(in   ) ::   kday     ! day value
1218      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldstop   ! stop if open to read a non-existing file (default = .TRUE.)
1219      !
1220      LOGICAL :: llprevyr              ! are we reading previous year  file?
1221      LOGICAL :: llprevmth             ! are we reading previous month file?
1222      INTEGER :: iyear, imonth, iday   ! first day of the current file in yyyy mm dd
1223      INTEGER :: isec_week             ! number of seconds since start of the weekly file
1224      INTEGER :: indexyr               ! year undex (O/1/2: previous/current/next)
1225      INTEGER :: iyear_len, imonth_len ! length (days) of iyear and imonth             !
1226      CHARACTER(len = 256)::   clname  ! temporary file name
1227      !!----------------------------------------------------------------------
1228      IF( PRESENT(kyear) ) THEN                             ! use given values
1229         iyear = kyear
1230         imonth = kmonth
1231         iday = kday
1232         IF ( sdjf%cltype(1:4) == 'week' ) THEN             ! find the day of the beginning of the week
1233            isec_week = ksec_week( sdjf%cltype(6:8) )- (86400 * 8 ) 
1234            llprevmth  = isec_week > nsec_month             ! longer time since beginning of the week than the month
1235            llprevyr   = llprevmth .AND. nmonth == 1
1236            iyear  = nyear  - COUNT((/llprevyr /))
1237            imonth = nmonth - COUNT((/llprevmth/)) + 12 * COUNT((/llprevyr /))
1238            iday   = nday   + nmonth_len(nmonth-1) * COUNT((/llprevmth/)) - isec_week / NINT(rday)
1239         ENDIF
1240      ELSE                                                  ! use current day values
1241         IF ( sdjf%cltype(1:4) == 'week' ) THEN             ! find the day of the beginning of the week
1242            isec_week  = ksec_week( sdjf%cltype(6:8) )      ! second since the beginning of the week
1243            llprevmth  = isec_week > nsec_month             ! longer time since beginning of the week than the month
1244            llprevyr   = llprevmth .AND. nmonth == 1
1245         ELSE
1246            isec_week  = 0
1247            llprevmth  = .FALSE.
1248            llprevyr   = .FALSE.
1249         ENDIF
1250         iyear  = nyear  - COUNT((/llprevyr /))
1251         imonth = nmonth - COUNT((/llprevmth/)) + 12 * COUNT((/llprevyr /))
1252         iday   = nday   + nmonth_len(nmonth-1) * COUNT((/llprevmth/)) - isec_week / NINT(rday)
1253      ENDIF
1254
1255      ! build the new filename if not climatological data
1256      clname=TRIM(sdjf%clrootname)
1257      !
1258      ! note that sdjf%ln_clim is is only acting on the presence of the year in the file name
1259      IF( .NOT. sdjf%ln_clim ) THEN   
1260                                         WRITE(clname, '(a,"_y",i4.4)' ) TRIM( sdjf%clrootname ), iyear    ! add year
1261         IF( sdjf%cltype /= 'yearly' )   WRITE(clname, '(a,"m" ,i2.2)' ) TRIM( clname          ), imonth   ! add month
1262      ELSE
1263         ! build the new filename if climatological data
1264         IF( sdjf%cltype /= 'yearly' )   WRITE(clname, '(a,"_m",i2.2)' ) TRIM( sdjf%clrootname ), imonth   ! add month
1265      ENDIF
1266      IF( sdjf%cltype == 'daily' .OR. sdjf%cltype(1:4) == 'week' ) &
1267            &                            WRITE(clname, '(a,"d" ,i2.2)' ) TRIM( clname          ), iday     ! add day
1268      !
1269      IF( TRIM(clname) /= TRIM(sdjf%clname) .OR. sdjf%num == 0 ) THEN   ! new file to be open
1270         !
1271         sdjf%clname = TRIM(clname)
1272         IF( sdjf%num /= 0 )   CALL iom_close( sdjf%num )   ! close file if already open
1273         CALL iom_open( sdjf%clname, sdjf%num, ldstop = ldstop, ldiof =  LEN(TRIM(sdjf%wgtname)) > 0 )
1274         !
1275         ! find the last record to be read -> update sdjf%nreclast
1276         indexyr = iyear - nyear + 1
1277         iyear_len = nyear_len( indexyr )
1278         SELECT CASE ( indexyr )
1279         CASE ( 0 )   ;   imonth_len = 31   ! previous year -> imonth = 12
1280         CASE ( 1 )   ;   imonth_len = nmonth_len(imonth) 
1281         CASE ( 2 )   ;   imonth_len = 31   ! next     year -> imonth = 1
1282         END SELECT
1283         !
1284         ! last record to be read in the current file
1285         IF    ( sdjf%nfreqh == -12 ) THEN                 ;   sdjf%nreclast = 1    !  yearly mean
1286         ELSEIF( sdjf%nfreqh ==  -1 ) THEN                                          ! monthly mean
1287            IF(     sdjf%cltype      == 'monthly' ) THEN   ;   sdjf%nreclast = 1
1288            ELSE                                           ;   sdjf%nreclast = 12
1289            ENDIF
1290         ELSE                                                                       ! higher frequency mean (in hours)
1291            IF(     sdjf%cltype      == 'monthly' ) THEN   ;   sdjf%nreclast = NINT( 24 * imonth_len / sdjf%nfreqh )
1292            ELSEIF( sdjf%cltype(1:4) == 'week'    ) THEN   ;   sdjf%nreclast = NINT( 24 * 7          / sdjf%nfreqh )
1293            ELSEIF( sdjf%cltype      == 'daily'   ) THEN   ;   sdjf%nreclast = NINT( 24              / sdjf%nfreqh )
1294            ELSE                                           ;   sdjf%nreclast = NINT( 24 * iyear_len  / sdjf%nfreqh )
1295            ENDIF
1296         ENDIF
1297         !
1298      ENDIF
1299      !
1300   END SUBROUTINE fld_clopn
1301
1302
1303   SUBROUTINE fld_fill( sdf, sdf_n, cdir, cdcaller, cdtitle, cdnam, knoprint )
1304      !!---------------------------------------------------------------------
1305      !!                    ***  ROUTINE fld_fill  ***
1306      !!
1307      !! ** Purpose :   fill the data structure (sdf) with the associated information
1308      !!              read in namelist (sdf_n) and control print
1309      !!----------------------------------------------------------------------
1310      TYPE(FLD)  , DIMENSION(:)          , INTENT(inout) ::   sdf        ! structure of input fields (file informations, fields read)
1311      TYPE(FLD_N), DIMENSION(:)          , INTENT(in   ) ::   sdf_n      ! array of namelist information structures
1312      CHARACTER(len=*)                   , INTENT(in   ) ::   cdir       ! Root directory for location of flx files
1313      CHARACTER(len=*)                   , INTENT(in   ) ::   cdcaller   ! name of the calling routine
1314      CHARACTER(len=*)                   , INTENT(in   ) ::   cdtitle    ! description of the calling routine
1315      CHARACTER(len=*)                   , INTENT(in   ) ::   cdnam      ! name of the namelist from which sdf_n comes
1316      INTEGER                  , OPTIONAL, INTENT(in   ) ::   knoprint   ! no calling routine information printed
1317      !
1318      INTEGER  ::   jf   ! dummy indices
1319      !!---------------------------------------------------------------------
1320      !
1321      DO jf = 1, SIZE(sdf)
1322         sdf(jf)%clrootname = TRIM( cdir )//TRIM( sdf_n(jf)%clname )
1323         sdf(jf)%clname     = "not yet defined"
1324         sdf(jf)%nfreqh     = sdf_n(jf)%nfreqh
1325         sdf(jf)%clvar      = sdf_n(jf)%clvar
1326         sdf(jf)%ln_tint    = sdf_n(jf)%ln_tint
1327         sdf(jf)%ln_clim    = sdf_n(jf)%ln_clim
1328         sdf(jf)%cltype     = sdf_n(jf)%cltype
1329         sdf(jf)%num        = -1
1330         sdf(jf)%wgtname    = " "
1331         IF( LEN( TRIM(sdf_n(jf)%wname) ) > 0 )   sdf(jf)%wgtname = TRIM( cdir )//TRIM( sdf_n(jf)%wname )
1332         sdf(jf)%lsmname = " "
1333         IF( LEN( TRIM(sdf_n(jf)%lname) ) > 0 )   sdf(jf)%lsmname = TRIM( cdir )//TRIM( sdf_n(jf)%lname )
1334         sdf(jf)%vcomp      = sdf_n(jf)%vcomp
1335         sdf(jf)%rotn(:)    = .TRUE.   ! pretend to be rotated -> won't try to rotate data before the first call to fld_get
1336         IF( sdf(jf)%cltype(1:4) == 'week' .AND. nn_leapy == 0  )   &
1337            &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdf(jf)%clrootname)//') needs nn_leapy = 1')
1338         IF( sdf(jf)%cltype(1:4) == 'week' .AND. sdf(jf)%ln_clim )   &
1339            &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdf(jf)%clrootname)//') needs ln_clim = .FALSE.')
1340         sdf(jf)%nreclast = -1 ! Set to non zero default value to avoid errors, is updated to meaningful value during fld_clopn
1341      END DO
1342      !
1343      IF(lwp) THEN      ! control print
1344         WRITE(numout,*)
1345         IF( .NOT.PRESENT( knoprint) ) THEN
1346            WRITE(numout,*) TRIM( cdcaller )//' : '//TRIM( cdtitle )
1347            WRITE(numout,*) (/ ('~', jf = 1, LEN_TRIM( cdcaller ) ) /)
1348         ENDIF
1349         WRITE(numout,*) '   fld_fill : fill data structure with information from namelist '//TRIM( cdnam )
1350         WRITE(numout,*) '   ~~~~~~~~'
1351         WRITE(numout,*) '      list of files and frequency (>0: in hours ; <0 in months)'
1352         DO jf = 1, SIZE(sdf)
1353            WRITE(numout,*) '      root filename: '  , TRIM( sdf(jf)%clrootname ), '   variable name: ', TRIM( sdf(jf)%clvar )
1354            WRITE(numout,*) '         frequency: '      ,       sdf(jf)%nfreqh      ,   &
1355               &                  '   time interp: '    ,       sdf(jf)%ln_tint     ,   &
1356               &                  '   climatology: '    ,       sdf(jf)%ln_clim
1357            WRITE(numout,*) '         weights: '        , TRIM( sdf(jf)%wgtname    ),   &
1358               &                  '   pairing: '        , TRIM( sdf(jf)%vcomp      ),   &
1359               &                  '   data type: '      ,       sdf(jf)%cltype      ,   &
1360               &                  '   land/sea mask:'   , TRIM( sdf(jf)%lsmname    )
1361            call flush(numout)
1362         END DO
1363      ENDIF
1364      !
1365   END SUBROUTINE fld_fill
1366
1367
1368   SUBROUTINE wgt_list( sd, kwgt, ldxios )
1369      !!---------------------------------------------------------------------
1370      !!                    ***  ROUTINE wgt_list  ***
1371      !!
1372      !! ** Purpose :   search array of WGTs and find a weights file entry,
1373      !!                or return a new one adding it to the end if new entry.
1374      !!                the weights data is read in and restructured (fld_weight)
1375      !!----------------------------------------------------------------------
1376      TYPE( FLD ), INTENT(in   ) ::   sd        ! field with name of weights file
1377      INTEGER    , INTENT(inout) ::   kwgt      ! index of weights
1378      LOGICAL    , INTENT(in), OPTIONAL :: ldxios
1379      !
1380      INTEGER ::   kw, nestid   ! local integer
1381      LOGICAL ::   found        ! local logical
1382      LOGICAL ::   luxios
1383      luxios = .FALSE.
1384      if(PRESENT(ldxios)) luxios = ldxios
1385      !!----------------------------------------------------------------------
1386      !
1387      !! search down linked list
1388      !! weights filename is either present or we hit the end of the list
1389      found = .FALSE.
1390      !
1391      !! because agrif nest part of filenames are now added in iom_open
1392      !! to distinguish between weights files on the different grids, need to track
1393      !! nest number explicitly
1394      nestid = 0
1395#if defined key_agrif
1396      nestid = Agrif_Fixed()
1397#endif
1398      DO kw = 1, nxt_wgt-1
1399         IF( TRIM(ref_wgts(kw)%wgtname) == TRIM(sd%wgtname) .AND. &
1400             ref_wgts(kw)%nestid == nestid) THEN
1401            kwgt = kw
1402            found = .TRUE.
1403            EXIT
1404         ENDIF
1405      END DO
1406      IF( .NOT.found ) THEN
1407         kwgt = nxt_wgt
1408         CALL fld_weight( sd, luxios )
1409      ENDIF
1410      !
1411   END SUBROUTINE wgt_list
1412
1413
1414   SUBROUTINE wgt_print( )
1415      !!---------------------------------------------------------------------
1416      !!                    ***  ROUTINE wgt_print  ***
1417      !!
1418      !! ** Purpose :   print the list of known weights
1419      !!----------------------------------------------------------------------
1420      INTEGER ::   kw   !
1421      !!----------------------------------------------------------------------
1422      !
1423      DO kw = 1, nxt_wgt-1
1424         WRITE(numout,*) 'weight file:  ',TRIM(ref_wgts(kw)%wgtname)
1425         WRITE(numout,*) '      ddims:  ',ref_wgts(kw)%ddims(1),ref_wgts(kw)%ddims(2)
1426         WRITE(numout,*) '     numwgt:  ',ref_wgts(kw)%numwgt
1427         WRITE(numout,*) '     jpiwgt:  ',ref_wgts(kw)%jpiwgt
1428         WRITE(numout,*) '     jpjwgt:  ',ref_wgts(kw)%jpjwgt
1429         WRITE(numout,*) '    botleft:  ',ref_wgts(kw)%botleft
1430         WRITE(numout,*) '   topright:  ',ref_wgts(kw)%topright
1431         IF( ref_wgts(kw)%cyclic ) THEN
1432            WRITE(numout,*) '       cyclical'
1433            IF( ref_wgts(kw)%overlap > 0 ) WRITE(numout,*) '              with overlap of ', ref_wgts(kw)%overlap
1434         ELSE
1435            WRITE(numout,*) '       not cyclical'
1436         ENDIF
1437         IF( ASSOCIATED(ref_wgts(kw)%data_wgt) )  WRITE(numout,*) '       allocated'
1438      END DO
1439      !
1440   END SUBROUTINE wgt_print
1441
1442
1443   SUBROUTINE fld_weight( sd, ldxios )
1444      !!---------------------------------------------------------------------
1445      !!                    ***  ROUTINE fld_weight  ***
1446      !!
1447      !! ** Purpose :   create a new WGT structure and fill in data from file,
1448      !!              restructuring as required
1449      !!----------------------------------------------------------------------
1450      TYPE( FLD ), INTENT(in) ::   sd   ! field with name of weights file
1451      LOGICAL, INTENT(in), OPTIONAL :: ldxios
1452      !!
1453      INTEGER ::   jn         ! dummy loop indices
1454      INTEGER ::   inum       ! local logical unit
1455      INTEGER ::   id         ! local variable id
1456      INTEGER ::   ipk        ! local vertical dimension
1457      INTEGER ::   zwrap      ! local integer
1458      LOGICAL ::   cyclical   !
1459      CHARACTER (len=5) ::   aname   !
1460      INTEGER , DIMENSION(:), ALLOCATABLE ::   ddims
1461      INTEGER,  DIMENSION(jpi,jpj) ::   data_src
1462      REAL(wp), DIMENSION(jpi,jpj) ::   data_tmp
1463      LOGICAL                      ::   luxios
1464      !!----------------------------------------------------------------------
1465      luxios = .FALSE.
1466      IF(PRESENT(ldxios)) luxios = ldxios
1467      !
1468      IF( nxt_wgt > tot_wgts ) THEN
1469        CALL ctl_stop("fld_weight: weights array size exceeded, increase tot_wgts")
1470      ENDIF
1471      !
1472      !! new weights file entry, add in extra information
1473      !! a weights file represents a 2D grid of a certain shape, so we assume that the current
1474      !! input data file is representative of all other files to be opened and processed with the
1475      !! current weights file
1476
1477      !! open input data file (non-model grid)
1478      CALL iom_open( sd%clname, inum, ldiof =  LEN(TRIM(sd%wgtname)) > 0 )
1479
1480      !! get dimensions
1481      IF ( SIZE(sd%fnow, 3) > 1 ) THEN
1482         ALLOCATE( ddims(4) )
1483      ELSE
1484         ALLOCATE( ddims(3) )
1485      ENDIF
1486      id = iom_varid( inum, sd%clvar, ddims )
1487
1488      !! close it
1489      CALL iom_close( inum )
1490
1491      !! now open the weights file
1492
1493      CALL iom_open ( sd%wgtname, inum )   ! interpolation weights
1494      IF ( inum > 0 ) THEN
1495         IF(luxios) THEN
1496            CALL iom_swap(cbwxios_context)
1497         ENDIF
1498         !! determine whether we have an east-west cyclic grid
1499         !! from global attribute called "ew_wrap" in the weights file
1500         !! note that if not found, iom_getatt returns -999 and cyclic with no overlap is assumed
1501         !! since this is the most common forcing configuration
1502
1503         CALL iom_getatt(inum, 'ew_wrap', zwrap)
1504         IF( zwrap >= 0 ) THEN
1505            cyclical = .TRUE.
1506         ELSE IF( zwrap == -999 ) THEN
1507            cyclical = .TRUE.
1508            zwrap = 0
1509         ELSE
1510            cyclical = .FALSE.
1511         ENDIF
1512
1513         ref_wgts(nxt_wgt)%ddims(1) = ddims(1)
1514         ref_wgts(nxt_wgt)%ddims(2) = ddims(2)
1515         ref_wgts(nxt_wgt)%wgtname = sd%wgtname
1516         ref_wgts(nxt_wgt)%overlap = zwrap
1517         ref_wgts(nxt_wgt)%cyclic = cyclical
1518         ref_wgts(nxt_wgt)%nestid = 0
1519#if defined key_agrif
1520         ref_wgts(nxt_wgt)%nestid = Agrif_Fixed()
1521#endif
1522         !! weights file is stored as a set of weights (wgt01->wgt04 or wgt01->wgt16)
1523         !! for each weight wgtNN there is an integer array srcNN which gives the point in
1524         !! the input data grid which is to be multiplied by the weight
1525         !! they are both arrays on the model grid so the result of the multiplication is
1526         !! added into an output array on the model grid as a running sum
1527
1528         !! two possible cases: bilinear (4 weights) or bicubic (16 weights)
1529         id = iom_varid(inum, 'src05', ldstop=.FALSE.)
1530         IF( id <= 0) THEN
1531            ref_wgts(nxt_wgt)%numwgt = 4
1532         ELSE
1533            ref_wgts(nxt_wgt)%numwgt = 16
1534         ENDIF
1535
1536         ALLOCATE( ref_wgts(nxt_wgt)%data_jpi(jpi,jpj,4) )
1537         ALLOCATE( ref_wgts(nxt_wgt)%data_jpj(jpi,jpj,4) )
1538         ALLOCATE( ref_wgts(nxt_wgt)%data_wgt(jpi,jpj,ref_wgts(nxt_wgt)%numwgt) )
1539
1540         DO jn = 1,4
1541            aname = ' '
1542            WRITE(aname,'(a3,i2.2)') 'src',jn
1543            data_tmp(:,:) = 0
1544            IF(luxios) THEN
1545               CALL iom_get ( inum, jpdom_data, aname//TRIM(sd%wgtname), data_tmp(:,:), ldxios = luxios )
1546            ELSE
1547               CALL iom_get ( inum, jpdom_data, aname, data_tmp(:,:) )
1548            ENDIF
1549            data_src(:,:) = INT(data_tmp(:,:))
1550            ref_wgts(nxt_wgt)%data_jpj(:,:,jn) = 1 + (data_src(:,:)-1) / ref_wgts(nxt_wgt)%ddims(1)
1551            ref_wgts(nxt_wgt)%data_jpi(:,:,jn) = data_src(:,:) - ref_wgts(nxt_wgt)%ddims(1)*(ref_wgts(nxt_wgt)%data_jpj(:,:,jn)-1)
1552         END DO
1553
1554         DO jn = 1, ref_wgts(nxt_wgt)%numwgt
1555            aname = ' '
1556            WRITE(aname,'(a3,i2.2)') 'wgt',jn
1557            ref_wgts(nxt_wgt)%data_wgt(:,:,jn) = 0.0
1558            IF(luxios) THEN
1559               CALL iom_get ( inum, jpdom_data, aname//TRIM(sd%wgtname), ref_wgts(nxt_wgt)%data_wgt(:,:,jn), ldxios = luxios )
1560            ELSE
1561               CALL iom_get ( inum, jpdom_data, aname, ref_wgts(nxt_wgt)%data_wgt(:,:,jn) )
1562            ENDIF
1563         END DO
1564         CALL iom_close (inum)
1565 
1566         ! find min and max indices in grid
1567         ref_wgts(nxt_wgt)%botleft(1) = MINVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1568         ref_wgts(nxt_wgt)%botleft(2) = MINVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1569         ref_wgts(nxt_wgt)%topright(1) = MAXVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1570         ref_wgts(nxt_wgt)%topright(2) = MAXVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1571
1572         ! and therefore dimensions of the input box
1573         ref_wgts(nxt_wgt)%jpiwgt = ref_wgts(nxt_wgt)%topright(1) - ref_wgts(nxt_wgt)%botleft(1) + 1
1574         ref_wgts(nxt_wgt)%jpjwgt = ref_wgts(nxt_wgt)%topright(2) - ref_wgts(nxt_wgt)%botleft(2) + 1
1575
1576         ! shift indexing of source grid
1577         ref_wgts(nxt_wgt)%data_jpi(:,:,:) = ref_wgts(nxt_wgt)%data_jpi(:,:,:) - ref_wgts(nxt_wgt)%botleft(1) + 1
1578         ref_wgts(nxt_wgt)%data_jpj(:,:,:) = ref_wgts(nxt_wgt)%data_jpj(:,:,:) - ref_wgts(nxt_wgt)%botleft(2) + 1
1579
1580         ! create input grid, give it a halo to allow gradient calculations
1581         ! SA: +3 stencil is a patch to avoid out-of-bound computation in some configuration.
1582         ! a more robust solution will be given in next release
1583         ipk =  SIZE(sd%fnow, 3)
1584         ALLOCATE( ref_wgts(nxt_wgt)%fly_dta(ref_wgts(nxt_wgt)%jpiwgt+3, ref_wgts(nxt_wgt)%jpjwgt+3 ,ipk) )
1585         IF( ref_wgts(nxt_wgt)%cyclic ) ALLOCATE( ref_wgts(nxt_wgt)%col(1,ref_wgts(nxt_wgt)%jpjwgt+3,ipk) )
1586         !
1587         nxt_wgt = nxt_wgt + 1
1588         !
1589         IF(luxios) CALL iom_swap(cxios_context)
1590      ELSE
1591         CALL ctl_stop( '    fld_weight : unable to read the file ' )
1592      ENDIF
1593
1594      DEALLOCATE (ddims )
1595      !
1596   END SUBROUTINE fld_weight
1597
1598
1599   SUBROUTINE apply_seaoverland( clmaskfile, zfieldo, jpi1_lsm, jpi2_lsm, jpj1_lsm,   &
1600      &                          jpj2_lsm, itmpi, itmpj, itmpz, rec1_lsm, recn_lsm )
1601      !!---------------------------------------------------------------------
1602      !!                    ***  ROUTINE apply_seaoverland  ***
1603      !!
1604      !! ** Purpose :   avoid spurious fluxes in coastal or near-coastal areas
1605      !!                due to the wrong usage of "land" values from the coarse
1606      !!                atmospheric model when spatial interpolation is required
1607      !!      D. Delrosso INGV         
1608      !!----------------------------------------------------------------------
1609      INTEGER,                   INTENT(in   ) :: itmpi,itmpj,itmpz                    ! lengths
1610      INTEGER,                   INTENT(in   ) :: jpi1_lsm,jpi2_lsm,jpj1_lsm,jpj2_lsm  ! temporary indices
1611      INTEGER, DIMENSION(3),     INTENT(in   ) :: rec1_lsm,recn_lsm                    ! temporary arrays for start and length
1612      REAL(wp),DIMENSION (:,:,:),INTENT(inout) :: zfieldo                              ! input/output array for seaoverland application
1613      CHARACTER (len=100),       INTENT(in   ) :: clmaskfile                           ! land/sea mask file name
1614      !
1615      INTEGER :: inum,jni,jnj,jnz,jc   ! local indices
1616      REAL(wp),DIMENSION (:,:,:),ALLOCATABLE :: zslmec1             ! local array for land point detection
1617      REAL(wp),DIMENSION (:,:),  ALLOCATABLE :: zfieldn   ! array of forcing field with undeff for land points
1618      REAL(wp),DIMENSION (:,:),  ALLOCATABLE :: zfield    ! array of forcing field
1619      !!---------------------------------------------------------------------
1620      !
1621      ALLOCATE ( zslmec1(itmpi,itmpj,itmpz), zfieldn(itmpi,itmpj), zfield(itmpi,itmpj) )
1622      !
1623      ! Retrieve the land sea mask data
1624      CALL iom_open( clmaskfile, inum )
1625      SELECT CASE( SIZE(zfieldo(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),3) )
1626      CASE(1)
1627         CALL iom_get( inum, jpdom_unknown, 'LSM', zslmec1(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,1), 1, rec1_lsm, recn_lsm)
1628      CASE DEFAULT
1629         CALL iom_get( inum, jpdom_unknown, 'LSM', zslmec1(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:), 1, rec1_lsm, recn_lsm)
1630      END SELECT
1631      CALL iom_close( inum )
1632      !
1633      DO jnz=1,rec1_lsm(3)             !! Loop over k dimension
1634         !
1635         DO jni = 1, itmpi                               !! copy the original field into a tmp array
1636            DO jnj = 1, itmpj                            !! substituting undeff over land points
1637               zfieldn(jni,jnj) = zfieldo(jni,jnj,jnz)
1638               IF( zslmec1(jni,jnj,jnz) == 1. )   zfieldn(jni,jnj) = undeff_lsm
1639            END DO
1640         END DO
1641         !
1642         CALL seaoverland( zfieldn, itmpi, itmpj, zfield )
1643         DO jc = 1, nn_lsm
1644            CALL seaoverland( zfield, itmpi, itmpj, zfield )
1645         END DO
1646         !
1647         !   Check for Undeff and substitute original values
1648         IF( ANY(zfield==undeff_lsm) ) THEN
1649            DO jni = 1, itmpi
1650               DO jnj = 1, itmpj
1651                  IF( zfield(jni,jnj)==undeff_lsm )   zfield(jni,jnj) = zfieldo(jni,jnj,jnz)
1652               END DO
1653            END DO
1654         ENDIF
1655         !
1656         zfieldo(:,:,jnz) = zfield(:,:)
1657         !
1658      END DO                           !! End Loop over k dimension
1659      !
1660      DEALLOCATE ( zslmec1, zfieldn, zfield )
1661      !
1662   END SUBROUTINE apply_seaoverland 
1663
1664
1665   SUBROUTINE seaoverland( zfieldn, ileni, ilenj, zfield )
1666      !!---------------------------------------------------------------------
1667      !!                    ***  ROUTINE seaoverland  ***
1668      !!
1669      !! ** Purpose :   create shifted matrices for seaoverland application 
1670      !!      D. Delrosso INGV
1671      !!----------------------------------------------------------------------
1672      INTEGER                      , INTENT(in   ) :: ileni,ilenj   ! lengths
1673      REAL, DIMENSION (ileni,ilenj), INTENT(in   ) :: zfieldn       ! array of forcing field with undeff for land points
1674      REAL, DIMENSION (ileni,ilenj), INTENT(  out) :: zfield        ! array of forcing field
1675      !
1676      REAL   , DIMENSION (ileni,ilenj)   :: zmat1, zmat2, zmat3, zmat4  ! local arrays
1677      REAL   , DIMENSION (ileni,ilenj)   :: zmat5, zmat6, zmat7, zmat8  !   -     -
1678      REAL   , DIMENSION (ileni,ilenj)   :: zlsm2d                      !   -     -
1679      REAL   , DIMENSION (ileni,ilenj,8) :: zlsm3d                      !   -     -
1680      LOGICAL, DIMENSION (ileni,ilenj,8) :: ll_msknan3d                 ! logical mask for undeff detection
1681      LOGICAL, DIMENSION (ileni,ilenj)   :: ll_msknan2d                 ! logical mask for undeff detection
1682      !!----------------------------------------------------------------------
1683      zmat8 = eoshift( zfieldn , SHIFT=-1 , BOUNDARY = (/zfieldn(:,1)/)     , DIM=2 )
1684      zmat1 = eoshift( zmat8   , SHIFT=-1 , BOUNDARY = (/zmat8(1,:)/)       , DIM=1 )
1685      zmat2 = eoshift( zfieldn , SHIFT=-1 , BOUNDARY = (/zfieldn(1,:)/)     , DIM=1 )
1686      zmat4 = eoshift( zfieldn , SHIFT= 1 , BOUNDARY = (/zfieldn(:,ilenj)/) , DIM=2 )
1687      zmat3 = eoshift( zmat4   , SHIFT=-1 , BOUNDARY = (/zmat4(1,:)/)       , DIM=1 )
1688      zmat5 = eoshift( zmat4   , SHIFT= 1 , BOUNDARY = (/zmat4(ileni,:)/)   , DIM=1 )
1689      zmat6 = eoshift( zfieldn , SHIFT= 1 , BOUNDARY = (/zfieldn(ileni,:)/) , DIM=1 )
1690      zmat7 = eoshift( zmat8   , SHIFT= 1 , BOUNDARY = (/zmat8(ileni,:)/)   , DIM=1 )
1691      !
1692      zlsm3d  = RESHAPE( (/ zmat1, zmat2, zmat3, zmat4, zmat5, zmat6, zmat7, zmat8 /), (/ ileni, ilenj, 8 /))
1693      ll_msknan3d = .NOT.( zlsm3d  == undeff_lsm )
1694      ll_msknan2d = .NOT.( zfieldn == undeff_lsm )  ! FALSE where is Undeff (land)
1695      zlsm2d = SUM( zlsm3d, 3 , ll_msknan3d ) / MAX( 1 , COUNT( ll_msknan3d , 3 ) )
1696      WHERE( COUNT( ll_msknan3d , 3 ) == 0._wp )   zlsm2d = undeff_lsm
1697      zfield = MERGE( zfieldn, zlsm2d, ll_msknan2d )
1698      !
1699   END SUBROUTINE seaoverland
1700
1701
1702   SUBROUTINE fld_interp( num, clvar, kw, kk, dta,  &
1703                          &         nrec, lsmfile)     
1704      !!---------------------------------------------------------------------
1705      !!                    ***  ROUTINE fld_interp  ***
1706      !!
1707      !! ** Purpose :   apply weights to input gridded data to create data
1708      !!                on model grid
1709      !!----------------------------------------------------------------------
1710      INTEGER                   , INTENT(in   ) ::   num     ! stream number
1711      CHARACTER(LEN=*)          , INTENT(in   ) ::   clvar   ! variable name
1712      INTEGER                   , INTENT(in   ) ::   kw      ! weights number
1713      INTEGER                   , INTENT(in   ) ::   kk      ! vertical dimension of kk
1714      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   dta     ! output field on model grid
1715      INTEGER                   , INTENT(in   ) ::   nrec    ! record number to read (ie time slice)
1716      CHARACTER(LEN=*)          , INTENT(in   ) ::   lsmfile ! land sea mask file name
1717      !
1718      INTEGER, DIMENSION(3) ::   rec1, recn           ! temporary arrays for start and length
1719      INTEGER, DIMENSION(3) ::   rec1_lsm, recn_lsm   ! temporary arrays for start and length in case of seaoverland
1720      INTEGER ::   ii_lsm1,ii_lsm2,ij_lsm1,ij_lsm2    ! temporary indices
1721      INTEGER ::   jk, jn, jm, jir, jjr               ! loop counters
1722      INTEGER ::   ni, nj                             ! lengths
1723      INTEGER ::   jpimin,jpiwid                      ! temporary indices
1724      INTEGER ::   jpimin_lsm,jpiwid_lsm              ! temporary indices
1725      INTEGER ::   jpjmin,jpjwid                      ! temporary indices
1726      INTEGER ::   jpjmin_lsm,jpjwid_lsm              ! temporary indices
1727      INTEGER ::   jpi1,jpi2,jpj1,jpj2                ! temporary indices
1728      INTEGER ::   jpi1_lsm,jpi2_lsm,jpj1_lsm,jpj2_lsm   ! temporary indices
1729      INTEGER ::   itmpi,itmpj,itmpz                     ! lengths
1730      REAL(wp),DIMENSION(:,:,:), ALLOCATABLE ::   ztmp_fly_dta                 ! local array of values on input grid     
1731      !!----------------------------------------------------------------------
1732      !
1733      !! for weighted interpolation we have weights at four corners of a box surrounding
1734      !! a model grid point, each weight is multiplied by a grid value (bilinear case)
1735      !! or by a grid value and gradients at the corner point (bicubic case)
1736      !! so we need to have a 4 by 4 subgrid surrounding each model point to cover both cases
1737
1738      !! sub grid from non-model input grid which encloses all grid points in this nemo process
1739      jpimin = ref_wgts(kw)%botleft(1)
1740      jpjmin = ref_wgts(kw)%botleft(2)
1741      jpiwid = ref_wgts(kw)%jpiwgt
1742      jpjwid = ref_wgts(kw)%jpjwgt
1743
1744      !! when reading in, expand this sub-grid by one halo point all the way round for calculating gradients
1745      rec1(1) = MAX( jpimin-1, 1 )
1746      rec1(2) = MAX( jpjmin-1, 1 )
1747      rec1(3) = 1
1748      recn(1) = MIN( jpiwid+2, ref_wgts(kw)%ddims(1)-rec1(1)+1 )
1749      recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
1750      recn(3) = kk
1751
1752      !! where we need to put it in the non-nemo grid fly_dta
1753      !! note that jpi1 and jpj1 only differ from 1 when jpimin and jpjmin are 1
1754      !! (ie at the extreme west or south of the whole input grid) and similarly for jpi2 and jpj2
1755      jpi1 = 2 + rec1(1) - jpimin
1756      jpj1 = 2 + rec1(2) - jpjmin
1757      jpi2 = jpi1 + recn(1) - 1
1758      jpj2 = jpj1 + recn(2) - 1
1759
1760
1761      IF( LEN( TRIM(lsmfile) ) > 0 ) THEN
1762      !! indeces for ztmp_fly_dta
1763      ! --------------------------
1764         rec1_lsm(1)=MAX(rec1(1)-nn_lsm,1)  ! starting index for enlarged external data, x direction
1765         rec1_lsm(2)=MAX(rec1(2)-nn_lsm,1)  ! starting index for enlarged external data, y direction
1766         rec1_lsm(3) = 1                    ! vertical dimension
1767         recn_lsm(1)=MIN(rec1(1)-rec1_lsm(1)+recn(1)+nn_lsm,ref_wgts(kw)%ddims(1)-rec1_lsm(1)) ! n points in x direction
1768         recn_lsm(2)=MIN(rec1(2)-rec1_lsm(2)+recn(2)+nn_lsm,ref_wgts(kw)%ddims(2)-rec1_lsm(2)) ! n points in y direction
1769         recn_lsm(3) = kk                   ! number of vertical levels in the input file
1770
1771      !  Avoid out of bound
1772         jpimin_lsm = MAX( rec1_lsm(1)+1, 1 )
1773         jpjmin_lsm = MAX( rec1_lsm(2)+1, 1 )
1774         jpiwid_lsm = MIN( recn_lsm(1)-2,ref_wgts(kw)%ddims(1)-rec1(1)+1)
1775         jpjwid_lsm = MIN( recn_lsm(2)-2,ref_wgts(kw)%ddims(2)-rec1(2)+1)
1776
1777         jpi1_lsm = 2+rec1_lsm(1)-jpimin_lsm
1778         jpj1_lsm = 2+rec1_lsm(2)-jpjmin_lsm
1779         jpi2_lsm = jpi1_lsm + recn_lsm(1) - 1
1780         jpj2_lsm = jpj1_lsm + recn_lsm(2) - 1
1781
1782
1783         itmpi=jpi2_lsm-jpi1_lsm+1
1784         itmpj=jpj2_lsm-jpj1_lsm+1
1785         itmpz=kk
1786         ALLOCATE(ztmp_fly_dta(itmpi,itmpj,itmpz))
1787         ztmp_fly_dta(:,:,:) = 0.0
1788         SELECT CASE( SIZE(ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),3) )
1789         CASE(1)
1790              CALL iom_get( num, jpdom_unknown, clvar, ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,1),   &
1791                 &                                                                nrec, rec1_lsm, recn_lsm)
1792         CASE DEFAULT
1793              CALL iom_get( num, jpdom_unknown, clvar, ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),   &
1794                 &                                                                nrec, rec1_lsm, recn_lsm)
1795         END SELECT
1796         CALL apply_seaoverland(lsmfile,ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),                  &
1797                 &                                      jpi1_lsm,jpi2_lsm,jpj1_lsm,jpj2_lsm,                  &
1798                 &                                      itmpi,itmpj,itmpz,rec1_lsm,recn_lsm)
1799
1800
1801         ! Relative indeces for remapping
1802         ii_lsm1 = (rec1(1)-rec1_lsm(1))+1
1803         ii_lsm2 = (ii_lsm1+recn(1))-1
1804         ij_lsm1 = (rec1(2)-rec1_lsm(2))+1
1805         ij_lsm2 = (ij_lsm1+recn(2))-1
1806
1807         ref_wgts(kw)%fly_dta(:,:,:) = 0.0
1808         ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:) = ztmp_fly_dta(ii_lsm1:ii_lsm2,ij_lsm1:ij_lsm2,:)
1809         DEALLOCATE(ztmp_fly_dta)
1810
1811      ELSE
1812         
1813         ref_wgts(kw)%fly_dta(:,:,:) = 0.0
1814         SELECT CASE( SIZE(ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:),3) )
1815         CASE(1)
1816              CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,1), nrec, rec1, recn)
1817         CASE DEFAULT
1818              CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:), nrec, rec1, recn)
1819         END SELECT
1820      ENDIF
1821     
1822
1823      !! first four weights common to both bilinear and bicubic
1824      !! data_jpi, data_jpj have already been shifted to (1,1) corresponding to botleft
1825      !! note that we have to offset by 1 into fly_dta array because of halo
1826      dta(:,:,:) = 0.0
1827      DO jk = 1,4
1828        DO jn = 1, jpj
1829          DO jm = 1,jpi
1830            ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1831            nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1832            dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk) * ref_wgts(kw)%fly_dta(ni+1,nj+1,:)
1833          END DO
1834        END DO
1835      END DO
1836
1837      IF (ref_wgts(kw)%numwgt .EQ. 16) THEN
1838
1839        !! fix up halo points that we couldnt read from file
1840        IF( jpi1 == 2 ) THEN
1841           ref_wgts(kw)%fly_dta(jpi1-1,:,:) = ref_wgts(kw)%fly_dta(jpi1,:,:)
1842        ENDIF
1843        IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
1844           ref_wgts(kw)%fly_dta(jpi2+1,:,:) = ref_wgts(kw)%fly_dta(jpi2,:,:)
1845        ENDIF
1846        IF( jpj1 == 2 ) THEN
1847           ref_wgts(kw)%fly_dta(:,jpj1-1,:) = ref_wgts(kw)%fly_dta(:,jpj1,:)
1848        ENDIF
1849        IF( jpj2 + jpjmin - 1 == ref_wgts(kw)%ddims(2)+1 .AND. jpj2 .lt. jpjwid+2 ) THEN
1850           ref_wgts(kw)%fly_dta(:,jpj2+1,:) = 2.0*ref_wgts(kw)%fly_dta(:,jpj2,:) - ref_wgts(kw)%fly_dta(:,jpj2-1,:)
1851        ENDIF
1852
1853        !! if data grid is cyclic we can do better on east-west edges
1854        !! but have to allow for whether first and last columns are coincident
1855        IF( ref_wgts(kw)%cyclic ) THEN
1856           rec1(2) = MAX( jpjmin-1, 1 )
1857           recn(1) = 1
1858           recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
1859           jpj1 = 2 + rec1(2) - jpjmin
1860           jpj2 = jpj1 + recn(2) - 1
1861           IF( jpi1 == 2 ) THEN
1862              rec1(1) = ref_wgts(kw)%ddims(1) - ref_wgts(kw)%overlap
1863              SELECT CASE( SIZE( ref_wgts(kw)%col(:,jpj1:jpj2,:),3) )
1864              CASE(1)
1865                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,1), nrec, rec1, recn)
1866              CASE DEFAULT
1867                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, rec1, recn)
1868              END SELECT     
1869              ref_wgts(kw)%fly_dta(jpi1-1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
1870           ENDIF
1871           IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
1872              rec1(1) = 1 + ref_wgts(kw)%overlap
1873              SELECT CASE( SIZE( ref_wgts(kw)%col(:,jpj1:jpj2,:),3) )
1874              CASE(1)
1875                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,1), nrec, rec1, recn)
1876              CASE DEFAULT
1877                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, rec1, recn)
1878              END SELECT
1879              ref_wgts(kw)%fly_dta(jpi2+1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
1880           ENDIF
1881        ENDIF
1882
1883        ! gradient in the i direction
1884        DO jk = 1,4
1885          DO jn = 1, jpj
1886            DO jm = 1,jpi
1887              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1888              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1889              dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+4) * 0.5 *         &
1890                               (ref_wgts(kw)%fly_dta(ni+2,nj+1,:) - ref_wgts(kw)%fly_dta(ni,nj+1,:))
1891            END DO
1892          END DO
1893        END DO
1894
1895        ! gradient in the j direction
1896        DO jk = 1,4
1897          DO jn = 1, jpj
1898            DO jm = 1,jpi
1899              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1900              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1901              dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+8) * 0.5 *         &
1902                               (ref_wgts(kw)%fly_dta(ni+1,nj+2,:) - ref_wgts(kw)%fly_dta(ni+1,nj,:))
1903            END DO
1904          END DO
1905        END DO
1906
1907         ! gradient in the ij direction
1908         DO jk = 1,4
1909            DO jn = 1, jpj
1910               DO jm = 1,jpi
1911                  ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1912                  nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1913                  dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+12) * 0.25 * ( &
1914                               (ref_wgts(kw)%fly_dta(ni+2,nj+2,:) - ref_wgts(kw)%fly_dta(ni  ,nj+2,:)) -   &
1915                               (ref_wgts(kw)%fly_dta(ni+2,nj  ,:) - ref_wgts(kw)%fly_dta(ni  ,nj  ,:)))
1916               END DO
1917            END DO
1918         END DO
1919         !
1920      END IF
1921      !
1922   END SUBROUTINE fld_interp
1923
1924
1925   FUNCTION ksec_week( cdday )
1926      !!---------------------------------------------------------------------
1927      !!                    ***  FUNCTION kshift_week ***
1928      !!
1929      !! ** Purpose :   return the first 3 letters of the first day of the weekly file
1930      !!---------------------------------------------------------------------
1931      CHARACTER(len=*), INTENT(in)   ::   cdday   ! first 3 letters of the first day of the weekly file
1932      !!
1933      INTEGER                        ::   ksec_week      ! output variable
1934      INTEGER                        ::   ijul, ishift   ! local integer
1935      CHARACTER(len=3),DIMENSION(7)  ::   cl_week 
1936      !!----------------------------------------------------------------------
1937      cl_week = (/"sun","sat","fri","thu","wed","tue","mon"/)
1938      DO ijul = 1, 7
1939         IF( cl_week(ijul) == TRIM(cdday) ) EXIT
1940      END DO
1941      IF( ijul .GT. 7 )   CALL ctl_stop( 'ksec_week: wrong day for sdjf%cltype(6:8): '//TRIM(cdday) )
1942      !
1943      ishift = ijul * NINT(rday)
1944      !
1945      ksec_week = nsec_week + ishift
1946      ksec_week = MOD( ksec_week, 7*NINT(rday) )
1947      !
1948   END FUNCTION ksec_week
1949
1950   !!======================================================================
1951END MODULE fldread
Note: See TracBrowser for help on using the repository browser.