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 branches/2012/dev_r3438_LOCEAN15_PISLOB/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2012/dev_r3438_LOCEAN15_PISLOB/NEMOGCM/NEMO/OPA_SRC/SBC/fldread.F90 @ 3556

Last change on this file since 3556 was 3556, checked in by cetlod, 11 years ago

branch:2012/dev_r3438_LOCEAN15_PISLOB add the Tropical Cyclones module

  • Property svn:keywords set to Id
File size: 67.9 KB
Line 
1MODULE fldread
2   !!======================================================================
3   !!                       ***  MODULE  fldread  ***
4   !! Ocean forcing:  read input field for surface boundary condition
5   !!=====================================================================
6   !! History :  2.0  !  06-2006  (S. Masson, G. Madec) Original code
7   !!                 !  05-2008  (S. Alderson) Modified for Interpolation in memory
8   !!                 !                         from input grid to model grid
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   fld_read      : read input fields used for the computation of the
13   !!                   surface boundary condition
14   !!----------------------------------------------------------------------
15   USE oce             ! ocean dynamics and tracers
16   USE dom_oce         ! ocean space and time domain
17   USE phycst          ! ???
18   USE in_out_manager  ! I/O manager
19   USE iom             ! I/O manager library
20   USE geo2ocean       ! for vector rotation on to model grid
21   USE lib_mpp         ! MPP library
22   USE wrk_nemo        ! work arrays
23   USE ioipsl, ONLY :   ymds2ju, ju2ymds   ! for calendar
24
25   IMPLICIT NONE
26   PRIVATE   
27 
28   PUBLIC   fld_map    ! routine called by tides_init
29
30   TYPE, PUBLIC ::   FLD_N      !: Namelist field informations
31      CHARACTER(len = 256) ::   clname      ! generic name of the NetCDF flux file
32      INTEGER              ::   nfreqh      ! frequency of each flux file
33      CHARACTER(len = 34)  ::   clvar       ! generic name of the variable in the NetCDF flux file
34      LOGICAL              ::   ln_tint     ! time interpolation or not (T/F)
35      LOGICAL              ::   ln_clim     ! climatology or not (T/F)
36      CHARACTER(len = 8)   ::   cltype      ! type of data file 'daily', 'monthly' or yearly'
37      CHARACTER(len = 34)  ::   wname       ! generic name of a NetCDF weights file to be used, blank if not
38      CHARACTER(len = 34)  ::   vcomp       ! symbolic component name if a vector that needs rotation
39      !                                     ! a string starting with "U" or "V" for each component   
40      !                                     ! chars 2 onwards identify which components go together 
41   END TYPE FLD_N
42
43   TYPE, PUBLIC ::   FLD        !: Input field related variables
44      CHARACTER(len = 256)            ::   clrootname   ! generic name of the NetCDF file
45      CHARACTER(len = 256)            ::   clname       ! current name of the NetCDF file
46      INTEGER                         ::   nfreqh       ! frequency of each flux file
47      CHARACTER(len = 34)             ::   clvar        ! generic name of the variable in the NetCDF flux file
48      LOGICAL                         ::   ln_tint      ! time interpolation or not (T/F)
49      LOGICAL                         ::   ln_clim      ! climatology or not (T/F)
50      CHARACTER(len = 8)              ::   cltype       ! type of data file 'daily', 'monthly' or yearly'
51      INTEGER                         ::   num          ! iom id of the jpfld files to be read
52      INTEGER , DIMENSION(2)          ::   nrec_b       ! before record (1: index, 2: second since Jan. 1st 00h of nit000 year)
53      INTEGER , DIMENSION(2)          ::   nrec_a       ! after  record (1: index, 2: second since Jan. 1st 00h of nit000 year)
54      REAL(wp) , ALLOCATABLE, DIMENSION(:,:,:  ) ::   fnow       ! input fields interpolated to now time step
55      REAL(wp) , ALLOCATABLE, DIMENSION(:,:,:,:) ::   fdta       ! 2 consecutive record of input fields
56      CHARACTER(len = 256)            ::   wgtname      ! current name of the NetCDF weight file acting as a key
57      !                                                 ! into the WGTLIST structure
58      CHARACTER(len = 34)             ::   vcomp        ! symbolic name for a vector component that needs rotation
59      LOGICAL                         ::   rotn         ! flag to indicate whether field has been rotated
60   END TYPE FLD
61
62   TYPE, PUBLIC ::   MAP_POINTER      !: Array of integer pointers to 1D arrays
63      INTEGER, POINTER   ::  ptr(:)
64   END TYPE MAP_POINTER
65
66!$AGRIF_DO_NOT_TREAT
67
68   !! keep list of all weights variables so they're only read in once
69   !! need to add AGRIF directives not to process this structure
70   !! also need to force wgtname to include AGRIF nest number
71   TYPE         ::   WGT        !: Input weights related variables
72      CHARACTER(len = 256)                    ::   wgtname      ! current name of the NetCDF weight file
73      INTEGER , DIMENSION(2)                  ::   ddims        ! shape of input grid
74      INTEGER , DIMENSION(2)                  ::   botleft      ! top left corner of box in input grid containing
75      !                                                         ! current processor grid
76      INTEGER , DIMENSION(2)                  ::   topright     ! top right corner of box
77      INTEGER                                 ::   jpiwgt       ! width of box on input grid
78      INTEGER                                 ::   jpjwgt       ! height of box on input grid
79      INTEGER                                 ::   numwgt       ! number of weights (4=bilinear, 16=bicubic)
80      INTEGER                                 ::   nestid       ! for agrif, keep track of nest we're in
81      INTEGER                                 ::   overlap      ! =0 when cyclic grid has no overlapping EW columns
82      !                                                         ! =>1 when they have one or more overlapping columns     
83      !                                                         ! =-1 not cyclic
84      LOGICAL                                 ::   cyclic       ! east-west cyclic or not
85      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpi     ! array of source integers
86      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpj     ! array of source integers
87      REAL(wp), DIMENSION(:,:,:), POINTER     ::   data_wgt     ! array of weights on model grid
88      REAL(wp), DIMENSION(:,:,:), POINTER     ::   fly_dta      ! array of values on input grid
89      REAL(wp), DIMENSION(:,:,:), POINTER     ::   col          ! temporary array for reading in columns
90   END TYPE WGT
91
92   INTEGER,     PARAMETER             ::   tot_wgts = 10
93   TYPE( WGT ), DIMENSION(tot_wgts)   ::   ref_wgts     ! array of wgts
94   INTEGER                            ::   nxt_wgt = 1  ! point to next available space in ref_wgts array
95
96!$AGRIF_END_DO_NOT_TREAT
97
98   PUBLIC   fld_read, fld_fill   ! called by sbc... modules
99
100   !!----------------------------------------------------------------------
101   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
102   !! $Id$
103   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
104   !!----------------------------------------------------------------------
105CONTAINS
106
107   SUBROUTINE fld_read( kt, kn_fsbc, sd, map, jit, time_offset )
108      !!---------------------------------------------------------------------
109      !!                    ***  ROUTINE fld_read  ***
110      !!                   
111      !! ** Purpose :   provide at each time step the surface ocean fluxes
112      !!                (momentum, heat, freshwater and runoff)
113      !!
114      !! ** Method  :   READ each input fields in NetCDF files using IOM
115      !!      and intepolate it to the model time-step.
116      !!         Several assumptions are made on the input file:
117      !!      blahblahblah....
118      !!----------------------------------------------------------------------
119      INTEGER  , INTENT(in   )               ::   kt        ! ocean time step
120      INTEGER  , INTENT(in   )               ::   kn_fsbc   ! sbc computation period (in time step)
121      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
122      TYPE(MAP_POINTER),INTENT(in), OPTIONAL, DIMENSION(:) ::   map   ! global-to-local mapping index
123      INTEGER  , INTENT(in   ), OPTIONAL     ::   jit       ! subcycle timestep for timesplitting option
124      INTEGER  , INTENT(in   ), OPTIONAL     ::   time_offset ! provide fields at time other than "now"
125                                                              ! time_offset = -1 => fields at "before" time level
126                                                              ! time_offset = +1 => fields at "after" time levels
127                                                              ! etc.
128      !!
129      INTEGER  ::   imf        ! size of the structure sd
130      INTEGER  ::   jf         ! dummy indices
131      INTEGER  ::   ireclast   ! last record to be read in the current year file
132      INTEGER  ::   isecend    ! number of second since Jan. 1st 00h of nit000 year at nitend
133      INTEGER  ::   isecsbc    ! number of seconds between Jan. 1st 00h of nit000 year and the middle of sbc time step
134      INTEGER  ::   itime_add  ! local time offset variable
135      LOGICAL  ::   llnxtyr    ! open next year  file?
136      LOGICAL  ::   llnxtmth   ! open next month file?
137      LOGICAL  ::   llstop     ! stop is the file does not exist
138      LOGICAL  ::   ll_firstcall ! true if this is the first call to fld_read for this set of fields
139      REAL(wp) ::   ztinta     ! ratio applied to after  records when doing time interpolation
140      REAL(wp) ::   ztintb     ! ratio applied to before records when doing time interpolation
141      CHARACTER(LEN=1000) ::   clfmt   ! write format
142      !!---------------------------------------------------------------------
143      ll_firstcall = .false.
144      IF( PRESENT(jit) ) THEN
145         IF(kt == nit000 .and. jit == 1) ll_firstcall = .true.
146      ELSE
147         IF(kt == nit000) ll_firstcall = .true.
148      ENDIF
149
150      itime_add = 0
151      IF( PRESENT(time_offset) ) itime_add = time_offset
152         
153      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
154      IF( present(jit) ) THEN 
155         ! ignore kn_fsbc in this case
156         isecsbc = nsec_year + nsec1jan000 + (jit+itime_add)*rdt/REAL(nn_baro,wp) 
157      ELSE
158         isecsbc = nsec_year + nsec1jan000 + NINT(0.5 * REAL(kn_fsbc - 1,wp) * rdttra(1)) + itime_add * rdttra(1)  ! middle of sbc time step
159      ENDIF
160      imf = SIZE( sd )
161      !
162      IF( ll_firstcall ) THEN                      ! initialization
163         IF( PRESENT(map) ) THEN
164            DO jf = 1, imf 
165               CALL fld_init( kn_fsbc, sd(jf), map(jf)%ptr )  ! read each before field (put them in after as they will be swapped)
166            END DO
167         ELSE
168            DO jf = 1, imf 
169               CALL fld_init( kn_fsbc, sd(jf) )       ! read each before field (put them in after as they will be swapped)
170            END DO
171         ENDIF
172         IF( lwp ) CALL wgt_print()                ! control print
173         CALL fld_rot( kt, sd )                    ! rotate vector fiels if needed
174      ENDIF
175      !                                            ! ====================================== !
176      IF( MOD( kt-1, kn_fsbc ) == 0 ) THEN         ! update field at each kn_fsbc time-step !
177         !                                         ! ====================================== !
178         !
179         DO jf = 1, imf                            ! ---   loop over field   --- !
180           
181            IF( isecsbc > sd(jf)%nrec_a(2) .OR. ll_firstcall ) THEN  ! read/update the after data?
182
183               IF( sd(jf)%ln_tint ) THEN                             ! swap before record field and informations
184                  sd(jf)%nrec_b(:) = sd(jf)%nrec_a(:)
185!CDIR COLLAPSE
186                  sd(jf)%fdta(:,:,:,1) = sd(jf)%fdta(:,:,:,2)
187               ENDIF
188
189               IF( PRESENT(jit) ) THEN
190                  CALL fld_rec( kn_fsbc, sd(jf), time_offset=itime_add, jit=jit )              ! update record informations
191               ELSE
192                  CALL fld_rec( kn_fsbc, sd(jf), time_offset=itime_add )                       ! update record informations
193               ENDIF
194
195               ! do we have to change the year/month/week/day of the forcing field??
196               IF( sd(jf)%ln_tint ) THEN
197                  ! if we do time interpolation we will need to open next year/month/week/day file before the end of the current
198                  ! 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)
199                  ! will be larger than the record number that should be read for current year/month/week/day
200
201                  ! last record to be read in the current file
202                  IF    ( sd(jf)%nfreqh == -12 ) THEN                 ;   ireclast = 1    !  yearly mean
203                  ELSEIF( sd(jf)%nfreqh ==  -1 ) THEN                                     ! monthly mean
204                     IF(     sd(jf)%cltype      == 'monthly' ) THEN   ;   ireclast = 1
205                     ELSE                                             ;   ireclast = 12
206                     ENDIF
207                  ELSE                                                                    ! higher frequency mean (in hours)
208                     IF(     sd(jf)%cltype      == 'monthly' ) THEN   ;   ireclast = 24 * nmonth_len(nmonth) / sd(jf)%nfreqh 
209                     ELSEIF( sd(jf)%cltype(1:4) == 'week'    ) THEN   ;   ireclast = 24 * 7                  / sd(jf)%nfreqh
210                     ELSEIF( sd(jf)%cltype      == 'daily'   ) THEN   ;   ireclast = 24                      / sd(jf)%nfreqh
211                     ELSE                                             ;   ireclast = 24 * nyear_len(     1 ) / sd(jf)%nfreqh 
212                     ENDIF
213                  ENDIF
214
215                  ! do we need next file data?
216                  IF( sd(jf)%nrec_a(1) > ireclast ) THEN
217
218                     sd(jf)%nrec_a(1) = 1              ! force to read the first record of the next file
219
220                     IF( .NOT. sd(jf)%ln_clim ) THEN   ! close the current file and open a new one.
221
222                        llnxtmth = sd(jf)%cltype == 'monthly' .OR. nday == nmonth_len(nmonth)      ! open next month file?
223                        llnxtyr  = sd(jf)%cltype == 'yearly'  .OR. (nmonth == 12 .AND. llnxtmth)   ! open next year  file?
224
225                        ! if the run finishes at the end of the current year/month/week/day, we will allow next
226                        ! year/month/week/day file to be not present. If the run continue further than the current
227                        ! year/month/week/day, next year/month/week/day file must exist
228                        isecend = nsec_year + nsec1jan000 + (nitend - kt) * NINT(rdttra(1))   ! second at the end of the run
229                        llstop = isecend > sd(jf)%nrec_a(2)                                   ! read more than 1 record of next year
230
231                        CALL fld_clopn( sd(jf), nyear  + COUNT((/llnxtyr /))                                           ,         &
232                           &                    nmonth + COUNT((/llnxtmth/)) - 12                 * COUNT((/llnxtyr /)),         &
233                           &                    nday   + 1                   - nmonth_len(nmonth) * COUNT((/llnxtmth/)), llstop )
234
235                        IF( sd(jf)%num <= 0 .AND. .NOT. llstop ) THEN    ! next year file does not exist
236                           CALL ctl_warn('next year/month/week/day file: '//TRIM(sd(jf)%clname)//     &
237                              &     ' not present -> back to current year/month/day')
238                           CALL fld_clopn( sd(jf), nyear, nmonth, nday )       ! back to the current year/month/day
239                           sd(jf)%nrec_a(1) = ireclast     ! force to read the last record to be read in the current year file
240                        ENDIF
241
242                     ENDIF
243                  ENDIF
244
245               ELSE
246                  ! if we are not doing time interpolation, we must change the year/month/week/day of the file just after
247                  ! switching to the NEW year/month/week/day. If it is the case, we are at the beginning of the
248                  ! year/month/week/day when calling fld_rec so sd(jf)%nrec_a(1) = 1
249                  IF( sd(jf)%nrec_a(1) == 1 .AND. .NOT. ( sd(jf)%ln_clim .AND. sd(jf)%cltype == 'yearly' ) )   &
250                     &   CALL fld_clopn( sd(jf), nyear, nmonth, nday )
251               ENDIF
252
253               ! read after data
254               IF( PRESENT(map) ) THEN
255                  CALL fld_get( sd(jf), map(jf)%ptr )
256               ELSE
257                  CALL fld_get( sd(jf) )
258               ENDIF
259
260            ENDIF
261         END DO                                    ! --- end loop over field --- !
262
263         CALL fld_rot( kt, sd )                    ! rotate vector fiels if needed
264
265         DO jf = 1, imf                            ! ---   loop over field   --- !
266            !
267            IF( sd(jf)%ln_tint ) THEN              ! temporal interpolation
268               IF(lwp .AND. kt - nit000 <= 100 ) THEN
269                  clfmt = "('fld_read: var ', a, ' kt = ', i8, ' (', f7.2,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
270                     &    "', records b/a: ', i4.4, '/', i4.4, ' (days ', f7.2,'/', f7.2, ')')"
271                  WRITE(numout, clfmt)  TRIM( sd(jf)%clvar ), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,   &           
272                     & 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
273                  WRITE(numout, *) 'itime_add is : ',itime_add
274               ENDIF
275               ! temporal interpolation weights
276               ztinta =  REAL( isecsbc - sd(jf)%nrec_b(2), wp ) / REAL( sd(jf)%nrec_a(2) - sd(jf)%nrec_b(2), wp )
277               ztintb =  1. - ztinta
278!CDIR COLLAPSE
279               sd(jf)%fnow(:,:,:) = ztintb * sd(jf)%fdta(:,:,:,1) + ztinta * sd(jf)%fdta(:,:,:,2)
280            ELSE   ! nothing to do...
281               IF(lwp .AND. kt - nit000 <= 100 ) THEN
282                  clfmt = "('fld_read: var ', a, ' kt = ', i8,' (', f7.2,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
283                     &    "', record: ', i4.4, ' (days ', f7.2, ' <-> ', f7.2, ')')"
284                  WRITE(numout, clfmt) TRIM(sd(jf)%clvar), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,    &
285                     &                 sd(jf)%nrec_a(1), REAL(sd(jf)%nrec_b(2),wp)/rday, REAL(sd(jf)%nrec_a(2),wp)/rday
286               ENDIF
287            ENDIF
288            !
289            IF( kt == nitend - kn_fsbc + 1 )   CALL iom_close( sd(jf)%num )   ! Close the input files
290
291         END DO                                    ! --- end loop over field --- !
292         !
293         !                                         ! ====================================== !
294      ENDIF                                        ! update field at each kn_fsbc time-step !
295      !                                            ! ====================================== !
296      !
297   END SUBROUTINE fld_read
298
299
300   SUBROUTINE fld_init( kn_fsbc, sdjf, map )
301      !!---------------------------------------------------------------------
302      !!                    ***  ROUTINE fld_init  ***
303      !!
304      !! ** Purpose :  - if time interpolation, read before data
305      !!               - open current year file
306      !!----------------------------------------------------------------------
307      INTEGER  , INTENT(in   ) ::   kn_fsbc   ! sbc computation period (in time step)
308      TYPE(FLD), INTENT(inout) ::   sdjf      ! input field related variables
309      INTEGER  , INTENT(in), OPTIONAL, DIMENSION(:) :: map ! global-to-local mapping indices
310      !!
311      LOGICAL :: llprevyr              ! are we reading previous year  file?
312      LOGICAL :: llprevmth             ! are we reading previous month file?
313      LOGICAL :: llprevweek            ! are we reading previous week  file?
314      LOGICAL :: llprevday             ! are we reading previous day   file?
315      LOGICAL :: llprev                ! llprevyr .OR. llprevmth .OR. llprevweek .OR. llprevday
316      INTEGER :: idvar                 ! variable id
317      INTEGER :: inrec                 ! number of record existing for this variable
318      INTEGER :: iyear, imonth, iday   ! first day of the current file in yyyy mm dd
319      INTEGER :: isec_week             ! number of seconds since start of the weekly file
320      CHARACTER(LEN=1000) ::   clfmt   ! write format
321      !!---------------------------------------------------------------------
322     
323      ! some default definitions...
324      sdjf%num = 0   ! default definition for non-opened file
325      IF( sdjf%ln_clim )   sdjf%clname = TRIM( sdjf%clrootname )   ! file name defaut definition, never change in this case
326      llprevyr   = .FALSE.
327      llprevmth  = .FALSE.
328      llprevweek = .FALSE.
329      llprevday  = .FALSE.
330      isec_week  = 0
331           
332      IF( sdjf%cltype(1:4) == 'week' .AND. nn_leapy == 0 )   &
333         &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdjf%clrootname)//') needs nn_leapy = 1')
334      IF( sdjf%cltype(1:4) == 'week' .AND. sdjf%ln_clim  )   &
335         &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdjf%clrootname)//') needs ln_clim = .FALSE.')
336
337      ! define record informations
338      CALL fld_rec( kn_fsbc, sdjf, ldbefore = .TRUE. )  ! return before values in sdjf%nrec_a (as we will swap it later)
339
340      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
341
342      IF( sdjf%ln_tint ) THEN ! we need to read the previous record and we will put it in the current record structure
343
344         IF( sdjf%nrec_a(1) == 0  ) THEN   ! we redefine record sdjf%nrec_a(1) with the last record of previous year file
345            IF    ( sdjf%nfreqh == -12 ) THEN   ! yearly mean
346               IF( sdjf%cltype == 'yearly' ) THEN             ! yearly file
347                  sdjf%nrec_a(1) = 1                                                       ! force to read the unique record
348                  llprevyr  = .NOT. sdjf%ln_clim                                           ! use previous year  file?
349               ELSE
350                  CALL ctl_stop( "fld_init: yearly mean file must be in a yearly type of file: "//TRIM(sdjf%clname) )
351               ENDIF
352            ELSEIF( sdjf%nfreqh ==  -1 ) THEN   ! monthly mean
353               IF( sdjf%cltype == 'monthly' ) THEN            ! monthly file
354                  sdjf%nrec_a(1) = 1                                                       ! force to read the unique record
355                  llprevmth = .TRUE.                                                       ! use previous month file?
356                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
357               ELSE                                           ! yearly file
358                  sdjf%nrec_a(1) = 12                                                      ! force to read december mean
359                  llprevyr = .NOT. sdjf%ln_clim                                            ! use previous year  file?
360               ENDIF
361            ELSE                                ! higher frequency mean (in hours)
362               IF    ( sdjf%cltype      == 'monthly' ) THEN   ! monthly file
363                  sdjf%nrec_a(1) = 24 * nmonth_len(nmonth-1) / sdjf%nfreqh                 ! last record of previous month
364                  llprevmth = .TRUE.                                                       ! use previous month file?
365                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
366               ELSEIF( sdjf%cltype(1:4) == 'week'    ) THEN   ! weekly file
367                  llprevweek = .TRUE.                                                      ! use previous week  file?
368                  sdjf%nrec_a(1) = 24 * 7 / sdjf%nfreqh                                    ! last record of previous week
369                  isec_week = NINT(rday) * 7                                               ! add a shift toward previous week
370               ELSEIF( sdjf%cltype      == 'daily'   ) THEN   ! daily file
371                  sdjf%nrec_a(1) = 24 / sdjf%nfreqh                                        ! last record of previous day
372                  llprevday = .TRUE.                                                       ! use previous day   file?
373                  llprevmth = llprevday .AND. nday   == 1                                  ! use previous month file?
374                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
375               ELSE                                           ! yearly file
376                  sdjf%nrec_a(1) = 24 * nyear_len(0) / sdjf%nfreqh                         ! last record of previous year
377                  llprevyr = .NOT. sdjf%ln_clim                                            ! use previous year  file?
378               ENDIF
379            ENDIF
380         ENDIF
381         IF ( sdjf%cltype(1:4) == 'week' ) THEN
382            isec_week = isec_week + ksec_week( sdjf%cltype(6:8) )   ! second since the beginning of the week
383            llprevmth = isec_week > nsec_month                      ! longer time since the beginning of the week than the month
384            llprevyr  = llprevmth .AND. nmonth == 1
385         ENDIF
386         llprev = llprevyr .OR. llprevmth .OR. llprevweek .OR. llprevday
387         !
388         iyear  = nyear  - COUNT((/llprevyr /))
389         imonth = nmonth - COUNT((/llprevmth/)) + 12 * COUNT((/llprevyr /))
390         iday   = nday   - COUNT((/llprevday/)) + nmonth_len(nmonth-1) * COUNT((/llprevmth/)) - isec_week / NINT(rday)
391         !
392         CALL fld_clopn( sdjf, iyear, imonth, iday, .NOT. llprev )
393
394         ! if previous year/month/day file does not exist, we switch to the current year/month/day
395         IF( llprev .AND. sdjf%num <= 0 ) THEN
396            CALL ctl_warn( 'previous year/month/week/day file: '//TRIM(sdjf%clname)//   &
397               &           ' not present -> back to current year/month/week/day' )
398            ! we force to read the first record of the current year/month/day instead of last record of previous year/month/day
399            llprev = .FALSE.
400            sdjf%nrec_a(1) = 1
401            CALL fld_clopn( sdjf, nyear, nmonth, nday )
402         ENDIF
403         
404         IF( llprev ) THEN   ! check if the last record sdjf%nrec_n(1) exists in the file
405            idvar = iom_varid( sdjf%num, sdjf%clvar )                                        ! id of the variable sdjf%clvar
406            IF( idvar <= 0 )   RETURN
407            inrec = iom_file( sdjf%num )%dimsz( iom_file( sdjf%num )%ndims(idvar), idvar )   ! size of the last dim of idvar
408            sdjf%nrec_a(1) = MIN( sdjf%nrec_a(1), inrec )   ! make sure we select an existing record
409         ENDIF
410
411         ! read before data
412         IF( PRESENT(map) ) THEN
413            CALL fld_get( sdjf, map )  ! read before values in after arrays(as we will swap it later)
414         ELSE
415            CALL fld_get( sdjf )  ! read before values in after arrays(as we will swap it later)
416         ENDIF
417
418         clfmt = "('fld_init : time-interpolation for ', a, ' read previous record = ', i4, ' at time = ', f7.2, ' days')"
419         IF(lwp) WRITE(numout, clfmt) TRIM(sdjf%clvar), sdjf%nrec_a(1), REAL(sdjf%nrec_a(2),wp)/rday
420
421         IF( llprev )   CALL iom_close( sdjf%num )          ! force to close previous year file (-> redefine sdjf%num to 0)
422
423      ENDIF
424
425      ! make sure current year/month/day file is opened
426      IF( sdjf%num <= 0 ) THEN
427         !
428         IF ( sdjf%cltype(1:4) == 'week' ) THEN
429            isec_week  = ksec_week( sdjf%cltype(6:8) )      ! second since the beginning of the week
430            llprevmth  = isec_week > nsec_month             ! longer time since beginning of the week than the month
431            llprevyr   = llprevmth .AND. nmonth == 1
432         ELSE
433            isec_week  = 0
434            llprevmth  = .FALSE.
435            llprevyr   = .FALSE.
436         ENDIF
437         !
438         iyear  = nyear  - COUNT((/llprevyr /))
439         imonth = nmonth - COUNT((/llprevmth/)) + 12 * COUNT((/llprevyr /))
440         iday   = nday   + nmonth_len(nmonth-1) * COUNT((/llprevmth/)) - isec_week / NINT(rday)
441         !
442         CALL fld_clopn( sdjf, iyear, imonth, iday )
443      ENDIF 
444      !
445   END SUBROUTINE fld_init
446
447
448   SUBROUTINE fld_rec( kn_fsbc, sdjf, ldbefore, jit, time_offset )
449      !!---------------------------------------------------------------------
450      !!                    ***  ROUTINE fld_rec  ***
451      !!
452      !! ** Purpose : Compute
453      !!              if sdjf%ln_tint = .TRUE.
454      !!                  nrec_a: record number and its time (nrec_b is obtained from nrec_a when swapping)
455      !!              if sdjf%ln_tint = .FALSE.
456      !!                  nrec_a(1): record number
457      !!                  nrec_b(2) and nrec_a(2): time of the beginning and end of the record (for print only)
458      !!----------------------------------------------------------------------
459      INTEGER  , INTENT(in   )           ::   kn_fsbc   ! sbc computation period (in time step)
460      TYPE(FLD), INTENT(inout)           ::   sdjf      ! input field related variables
461      LOGICAL  , INTENT(in   ), OPTIONAL ::   ldbefore  ! sent back before record values (default = .FALSE.)
462      INTEGER  , INTENT(in   ), OPTIONAL ::   jit       ! index of barotropic subcycle
463                                                        ! used only if sdjf%ln_tint = .TRUE.
464      INTEGER  , INTENT(in   ), OPTIONAL ::   time_offset  ! Offset of required time level compared to "now"
465                                                           ! time level in units of time steps.
466      !!
467      LOGICAL  ::   llbefore    ! local definition of ldbefore
468      INTEGER  ::   iendrec     ! end of this record (in seconds)
469      INTEGER  ::   imth        ! month number
470      INTEGER  ::   ifreq_sec   ! frequency mean (in seconds)
471      INTEGER  ::   isec_week   ! number of seconds since the start of the weekly file
472      INTEGER  ::   itime_add   ! local time offset variable
473      REAL(wp) ::   ztmp        ! temporary variable
474      !!----------------------------------------------------------------------
475      !
476      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
477      !
478      IF( PRESENT(ldbefore) ) THEN   ;   llbefore = ldbefore .AND. sdjf%ln_tint   ! needed only if sdjf%ln_tint = .TRUE.
479      ELSE                           ;   llbefore = .FALSE.
480      ENDIF
481      !
482      itime_add = 0
483      IF( PRESENT(time_offset) ) itime_add = time_offset
484      !
485      !                                      ! =========== !
486      IF    ( sdjf%nfreqh == -12 ) THEN      ! yearly mean
487         !                                   ! =========== !
488         !
489         IF( sdjf%ln_tint ) THEN                 ! time interpolation, shift by 1/2 record
490            !
491            !                  INT( ztmp )
492            !                     /|\
493            !                    1 |    *----
494            !                    0 |----(             
495            !                      |----+----|--> time
496            !                      0   /|\   1   (nday/nyear_len(1))
497            !                           |   
498            !                           |   
499            !       forcing record :    1
500            !                           
501            ztmp = REAL( nday, wp ) / REAL( nyear_len(1), wp ) + 0.5
502            IF( PRESENT(jit) ) THEN
503               ztmp = ztmp + (jit+itime_add)*rdt/REAL(nn_baro,wp)
504            ELSE
505               ztmp = ztmp + itime_add*rdttra(1)
506            ENDIF
507            sdjf%nrec_a(1) = 1 + INT( ztmp ) - COUNT((/llbefore/))
508            ! swap at the middle of the year
509            IF( llbefore ) THEN   ;   sdjf%nrec_a(2) = nsec1jan000 - NINT(0.5 * rday) * nyear_len(0)
510            ELSE                  ;   sdjf%nrec_a(2) = nsec1jan000 + NINT(0.5 * rday) * nyear_len(1)   
511            ENDIF
512         ELSE                                    ! no time interpolation
513            sdjf%nrec_a(1) = 1
514            sdjf%nrec_a(2) = NINT(rday) * nyear_len(1) + nsec1jan000   ! swap at the end    of the year
515            sdjf%nrec_b(2) = nsec1jan000                               ! beginning of the year (only for print)
516         ENDIF
517         !
518         !                                   ! ============ !
519      ELSEIF( sdjf%nfreqh ==  -1 ) THEN      ! monthly mean !
520         !                                   ! ============ !
521         !
522         IF( sdjf%ln_tint ) THEN                 ! time interpolation, shift by 1/2 record
523            !
524            !                  INT( ztmp )
525            !                     /|\
526            !                    1 |    *----
527            !                    0 |----(             
528            !                      |----+----|--> time
529            !                      0   /|\   1   (nday/nmonth_len(nmonth))
530            !                           |   
531            !                           |   
532            !       forcing record :  nmonth
533            !                           
534            ztmp = REAL( nday, wp ) / REAL( nmonth_len(nmonth), wp ) + 0.5
535            IF( PRESENT(jit) ) THEN
536               ztmp = ztmp + (jit+itime_add)*rdt/REAL(nn_baro,wp)
537            ELSE
538               ztmp = ztmp + itime_add*rdttra(1)
539            ENDIF
540            imth = nmonth + INT( ztmp ) - COUNT((/llbefore/))
541            IF( sdjf%cltype == 'monthly' ) THEN   ;   sdjf%nrec_a(1) = 1 + INT( ztmp ) - COUNT((/llbefore/))
542            ELSE                                  ;   sdjf%nrec_a(1) = imth
543            ENDIF
544            sdjf%nrec_a(2) = nmonth_half(   imth ) + nsec1jan000   ! swap at the middle of the month
545         ELSE                                    ! no time interpolation
546            IF( sdjf%cltype == 'monthly' ) THEN   ;   sdjf%nrec_a(1) = 1
547            ELSE                                  ;   sdjf%nrec_a(1) = nmonth
548            ENDIF
549            sdjf%nrec_a(2) =  nmonth_end(nmonth  ) + nsec1jan000   ! swap at the end    of the month
550            sdjf%nrec_b(2) =  nmonth_end(nmonth-1) + nsec1jan000   ! beginning of the month (only for print)
551         ENDIF
552         !
553         !                                   ! ================================ !
554      ELSE                                   ! higher frequency mean (in hours)
555         !                                   ! ================================ !
556         !
557         ifreq_sec = sdjf%nfreqh * 3600                                                 ! frequency mean (in seconds)
558         IF( sdjf%cltype(1:4) == 'week' )   isec_week = ksec_week( sdjf%cltype(6:8) )   ! since the first day of the current week
559         ! number of second since the beginning of the file
560         IF(     sdjf%cltype      == 'monthly' ) THEN   ;   ztmp = REAL(nsec_month,wp)  ! since the first day of the current month
561         ELSEIF( sdjf%cltype(1:4) == 'week'    ) THEN   ;   ztmp = REAL(isec_week ,wp)  ! since the first day of the current week
562         ELSEIF( sdjf%cltype      == 'daily'   ) THEN   ;   ztmp = REAL(nsec_day  ,wp)  ! since 00h of the current day
563         ELSE                                           ;   ztmp = REAL(nsec_year ,wp)  ! since 00h on Jan 1 of the current year
564         ENDIF
565         ztmp = ztmp + 0.5 * REAL(kn_fsbc - 1, wp) * rdttra(1)   ! shift time to be centrered in the middle of sbc time step
566         ztmp = ztmp + 0.01 * rdttra(1)                          ! add 0.01 time step to avoid truncation error
567         IF( PRESENT(jit) ) THEN
568            ztmp = ztmp + (jit+itime_add)*rdt/REAL(nn_baro,wp)
569         ELSE
570            ztmp = ztmp + itime_add*rdttra(1)
571         ENDIF
572         IF( sdjf%ln_tint ) THEN                ! time interpolation, shift by 1/2 record
573            !
574            !                  INT( ztmp )
575            !                     /|\
576            !                    2 |        *-----(
577            !                    1 |  *-----(
578            !                    0 |--(             
579            !                      |--+--|--+--|--+--|--> time
580            !                      0 /|\ 1 /|\ 2 /|\ 3 (nsec_year/ifreq_sec) or (nsec_month/ifreq_sec)
581            !                         |     |     |
582            !                         |     |     |
583            !       forcing record :  1     2     3
584            !                   
585            ztmp= ztmp / REAL(ifreq_sec, wp) + 0.5
586         ELSE                                   ! no time interpolation
587            !
588            !                  INT( ztmp )
589            !                     /|\
590            !                    2 |           *-----(
591            !                    1 |     *-----(
592            !                    0 |-----(             
593            !                      |--+--|--+--|--+--|--> time
594            !                      0 /|\ 1 /|\ 2 /|\ 3 (nsec_year/ifreq_sec) or (nsec_month/ifreq_sec)
595            !                         |     |     |
596            !                         |     |     |
597            !       forcing record :  1     2     3
598            !                           
599            ztmp= ztmp / REAL(ifreq_sec, wp)
600         ENDIF
601         sdjf%nrec_a(1) = 1 + INT( ztmp ) - COUNT((/llbefore/))   ! record nomber to be read
602
603         iendrec = ifreq_sec * sdjf%nrec_a(1) + nsec1jan000       ! end of this record (in second)
604         ! add the number of seconds between 00h Jan 1 and the end of previous month/week/day (ok if nmonth=1)
605         IF( sdjf%cltype      == 'monthly' )   iendrec = iendrec + NINT(rday) * SUM(nmonth_len(1:nmonth -1))
606         IF( sdjf%cltype(1:4) == 'week'    )   iendrec = iendrec + ( nsec_year - isec_week )
607         IF( sdjf%cltype      == 'daily'   )   iendrec = iendrec + NINT(rday) * ( nday_year - 1 )
608         IF( sdjf%ln_tint ) THEN
609             sdjf%nrec_a(2) = iendrec - ifreq_sec / 2        ! swap at the middle of the record
610         ELSE
611             sdjf%nrec_a(2) = iendrec                        ! swap at the end    of the record
612             sdjf%nrec_b(2) = iendrec - ifreq_sec            ! beginning of the record (only for print)
613         ENDIF
614         !
615      ENDIF
616      !
617   END SUBROUTINE fld_rec
618
619
620   SUBROUTINE fld_get( sdjf, map )
621      !!---------------------------------------------------------------------
622      !!                    ***  ROUTINE fld_get  ***
623      !!
624      !! ** Purpose :   read the data
625      !!----------------------------------------------------------------------
626      TYPE(FLD), INTENT(inout) ::   sdjf   ! input field related variables
627      INTEGER  , INTENT(in), OPTIONAL, DIMENSION(:) :: map ! global-to-local mapping indices
628      !!
629      INTEGER                  ::   ipi    ! number of i-point sdjf%fdta
630      INTEGER                  ::   ipk    ! number of vertical levels of sdjf%fdta ( 2D: ipk=1 ; 3D: ipk=jpk )
631      INTEGER                  ::   iw     ! index into wgts array
632      !!---------------------------------------------------------------------
633           
634      ipi = SIZE( sdjf%fnow, 1 )
635      ipk = SIZE( sdjf%fnow, 3 )
636
637      IF( PRESENT(map) ) THEN
638         IF( sdjf%ln_tint ) THEN   ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fdta(:,:,:,2), sdjf%nrec_a(1), map )
639         ELSE                      ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fnow(:,:,:  ), sdjf%nrec_a(1), map )
640         ENDIF
641      ELSE IF( LEN(TRIM(sdjf%wgtname)) > 0 ) THEN
642         CALL wgt_list( sdjf, iw )
643         IF( sdjf%ln_tint ) THEN   ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw , ipk  , sdjf%fdta(:,:,:,2), sdjf%nrec_a(1) )
644         ELSE                      ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw , ipk  , sdjf%fnow(:,:,:  ), sdjf%nrec_a(1) )
645         ENDIF
646      ELSE
647         SELECT CASE( ipk )
648         CASE(1)   
649            IF( ipi == jpi ) THEN
650               IF( sdjf%ln_tint ) THEN   ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fdta(:,:,1,2), sdjf%nrec_a(1) )
651               ELSE                      ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fnow(:,:,1  ), sdjf%nrec_a(1) )
652               ENDIF
653            ELSE
654               IF( sdjf%ln_tint ) THEN   ;   CALL iom_get( sdjf%num, jpdom_unknown, sdjf%clvar, sdjf%fdta(:,:,1,2), sdjf%nrec_a(1) )
655               ELSE                      ;   CALL iom_get( sdjf%num, jpdom_unknown, sdjf%clvar, sdjf%fnow(:,:,1  ), sdjf%nrec_a(1) )
656               ENDIF
657            ENDIF
658         CASE DEFAULT
659            IF( sdjf%ln_tint ) THEN   ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fdta(:,:,:,2), sdjf%nrec_a(1) )
660            ELSE                      ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fnow(:,:,:  ), sdjf%nrec_a(1) )
661            ENDIF
662         END SELECT
663      ENDIF
664      !
665      sdjf%rotn = .false.   ! vector not yet rotated
666
667   END SUBROUTINE fld_get
668
669   SUBROUTINE fld_map( num, clvar, dta, nrec, map )
670      !!---------------------------------------------------------------------
671      !!                    ***  ROUTINE fld_get  ***
672      !!
673      !! ** Purpose :   read global data from file and map onto local data
674      !!                using a general mapping (for open boundaries)
675      !!----------------------------------------------------------------------
676#if defined key_bdy
677      USE bdy_oce, ONLY:  dta_global         ! workspace to read in global data arrays
678#endif
679
680      INTEGER                   , INTENT(in ) ::   num     ! stream number
681      CHARACTER(LEN=*)          , INTENT(in ) ::   clvar   ! variable name
682      REAL(wp), DIMENSION(:,:,:), INTENT(out) ::   dta   ! output field on model grid (2 dimensional)
683      INTEGER                   , INTENT(in ) ::   nrec    ! record number to read (ie time slice)
684      INTEGER,  DIMENSION(:)    , INTENT(in ) ::   map     ! global-to-local mapping indices
685      !!
686      INTEGER                                 ::   ipi      ! length of boundary data on local process
687      INTEGER                                 ::   ipj      ! length of dummy dimension ( = 1 )
688      INTEGER                                 ::   ipk      ! number of vertical levels of dta ( 2D: ipk=1 ; 3D: ipk=jpk )
689      INTEGER                                 ::   ilendta  ! length of data in file
690      INTEGER                                 ::   idvar    ! variable ID
691      INTEGER                                 ::   ib, ik   ! loop counters
692      INTEGER                                 ::   ierr
693      REAL(wp), POINTER, DIMENSION(:,:,:)     ::   dta_read ! work space for global data
694      !!---------------------------------------------------------------------
695           
696#if defined key_bdy
697      dta_read => dta_global
698#endif
699
700      ipi = SIZE( dta, 1 )
701      ipj = 1
702      ipk = SIZE( dta, 3 )
703
704      idvar   = iom_varid( num, clvar )
705      ilendta = iom_file(num)%dimsz(1,idvar)
706      IF(lwp) WRITE(numout,*) 'Dim size for ',TRIM(clvar),' is ', ilendta
707      IF(lwp) WRITE(numout,*) 'Number of levels for ',TRIM(clvar),' is ', ipk
708
709      SELECT CASE( ipk )
710      CASE(1)   
711         CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1    ), nrec )
712      CASE DEFAULT
713         CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1:ipk), nrec )
714      END SELECT
715      !
716      DO ib = 1, ipi
717         DO ik = 1, ipk
718            dta(ib,1,ik) =  dta_read(map(ib),1,ik)
719         END DO
720      END DO
721
722   END SUBROUTINE fld_map
723
724
725   SUBROUTINE fld_rot( kt, sd )
726      !!---------------------------------------------------------------------
727      !!                    ***  ROUTINE fld_rot  ***
728      !!
729      !! ** Purpose :   Vector fields may need to be rotated onto the local grid direction
730      !!----------------------------------------------------------------------
731      INTEGER  , INTENT(in   )               ::   kt        ! ocean time step
732      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
733      !!
734      INTEGER                           ::   ju, jv, jk   ! loop indices
735      INTEGER                           ::   imf          ! size of the structure sd
736      INTEGER                           ::   ill          ! character length
737      INTEGER                           ::   iv           ! indice of V component
738      REAL(wp), POINTER, DIMENSION(:,:) ::   utmp, vtmp   ! temporary arrays for vector rotation
739      CHARACTER (LEN=100)               ::   clcomp       ! dummy weight name
740      !!---------------------------------------------------------------------
741
742      CALL wrk_alloc( jpi,jpj, utmp, vtmp )
743
744      !! (sga: following code should be modified so that pairs arent searched for each time
745      !
746      imf = SIZE( sd )
747      DO ju = 1, imf
748         ill = LEN_TRIM( sd(ju)%vcomp )
749         IF( ill > 0 .AND. .NOT. sd(ju)%rotn ) THEN   ! find vector rotations required             
750             IF( sd(ju)%vcomp(1:1) == 'U' ) THEN      ! east-west component has symbolic name starting with 'U'
751                ! look for the north-south component which has same symbolic name but with 'U' replaced with 'V'
752                clcomp = 'V' // sd(ju)%vcomp(2:ill)   ! works even if ill == 1
753                iv = -1
754                DO jv = 1, imf
755                  IF( TRIM(sd(jv)%vcomp) == TRIM(clcomp) )   iv = jv
756                END DO
757                IF( iv > 0 ) THEN   ! fields ju and iv are two components which need to be rotated together
758                   DO jk = 1, SIZE( sd(ju)%fnow, 3 )
759                      IF( sd(ju)%ln_tint )THEN
760                         CALL rot_rep( sd(ju)%fdta(:,:,jk,2), sd(iv)%fdta(:,:,jk,2), 'T', 'en->i', utmp(:,:) )
761                         CALL rot_rep( sd(ju)%fdta(:,:,jk,2), sd(iv)%fdta(:,:,jk,2), 'T', 'en->j', vtmp(:,:) )
762                         sd(ju)%fdta(:,:,jk,2) = utmp(:,:)   ;   sd(iv)%fdta(:,:,jk,2) = vtmp(:,:)
763                      ELSE
764                         CALL rot_rep( sd(ju)%fnow(:,:,jk  ), sd(iv)%fnow(:,:,jk  ), 'T', 'en->i', utmp(:,:) )
765                         CALL rot_rep( sd(ju)%fnow(:,:,jk  ), sd(iv)%fnow(:,:,jk  ), 'T', 'en->j', vtmp(:,:) )
766                         sd(ju)%fnow(:,:,jk  ) = utmp(:,:)   ;   sd(iv)%fnow(:,:,jk  ) = vtmp(:,:)
767                      ENDIF
768                   END DO
769                   sd(ju)%rotn = .TRUE.               ! vector was rotated
770                   IF( lwp .AND. kt == nit000 )   WRITE(numout,*)   &
771                      &   'fld_read: vector pair ('//TRIM(sd(ju)%clvar)//', '//TRIM(sd(iv)%clvar)//') rotated on to model grid'
772                ENDIF
773             ENDIF
774          ENDIF
775       END DO
776      !
777      CALL wrk_dealloc( jpi,jpj, utmp, vtmp )
778      !
779   END SUBROUTINE fld_rot
780
781
782   SUBROUTINE fld_clopn( sdjf, kyear, kmonth, kday, ldstop )
783      !!---------------------------------------------------------------------
784      !!                    ***  ROUTINE fld_clopn  ***
785      !!
786      !! ** Purpose :   update the file name and open the file
787      !!----------------------------------------------------------------------
788      TYPE(FLD)        , INTENT(inout) ::   sdjf     ! input field related variables
789      INTEGER          , INTENT(in   ) ::   kyear    ! year value
790      INTEGER          , INTENT(in   ) ::   kmonth   ! month value
791      INTEGER          , INTENT(in   ) ::   kday     ! day value
792      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldstop   ! stop if open to read a non-existing file (default = .TRUE.)
793      !!----------------------------------------------------------------------
794
795      IF( sdjf%num /= 0 )   CALL iom_close( sdjf%num )   ! close file if already open
796      ! build the new filename if not climatological data
797      sdjf%clname=TRIM(sdjf%clrootname)
798      !
799      ! note that sdjf%ln_clim is is only acting on presence of the year in the file
800      IF( .NOT. sdjf%ln_clim ) THEN   
801                                         WRITE(sdjf%clname, '(a,"_y",i4.4)' ) TRIM( sdjf%clrootname ), kyear    ! add year
802         IF( sdjf%cltype /= 'yearly' )   WRITE(sdjf%clname, '(a,"m" ,i2.2)' ) TRIM( sdjf%clname     ), kmonth   ! add month
803      ELSE
804         ! build the new filename if climatological data
805         IF( sdjf%cltype /= 'yearly' )   WRITE(sdjf%clname, '(a,"_m",i2.2)' ) TRIM( sdjf%clrootname ), kmonth   ! add month
806      ENDIF
807      IF( sdjf%cltype == 'daily' .OR. sdjf%cltype(1:4) == 'week' ) &
808            &                            WRITE(sdjf%clname, '(a,"d" ,i2.2)' ) TRIM( sdjf%clname     ), kday     ! add day
809      !
810      CALL iom_open( sdjf%clname, sdjf%num, ldstop = ldstop, ldiof =  LEN(TRIM(sdjf%wgtname)) > 0 )
811     !
812   END SUBROUTINE fld_clopn
813
814
815   SUBROUTINE fld_fill( sdf, sdf_n, cdir, cdcaller, cdtitle, cdnam )
816      !!---------------------------------------------------------------------
817      !!                    ***  ROUTINE fld_fill  ***
818      !!
819      !! ** Purpose :   fill sdf with sdf_n and control print
820      !!----------------------------------------------------------------------
821      TYPE(FLD)  , DIMENSION(:), INTENT(inout) ::   sdf        ! structure of input fields (file informations, fields read)
822      TYPE(FLD_N), DIMENSION(:), INTENT(in   ) ::   sdf_n      ! array of namelist information structures
823      CHARACTER(len=*)         , INTENT(in   ) ::   cdir       ! Root directory for location of flx files
824      CHARACTER(len=*)         , INTENT(in   ) ::   cdcaller   !
825      CHARACTER(len=*)         , INTENT(in   ) ::   cdtitle    !
826      CHARACTER(len=*)         , INTENT(in   ) ::   cdnam      !
827      !
828      INTEGER  ::   jf       ! dummy indices
829      !!---------------------------------------------------------------------
830
831      DO jf = 1, SIZE(sdf)
832         sdf(jf)%clrootname = TRIM( cdir )//TRIM( sdf_n(jf)%clname )
833         sdf(jf)%nfreqh     = sdf_n(jf)%nfreqh
834         sdf(jf)%clvar      = sdf_n(jf)%clvar
835         sdf(jf)%ln_tint    = sdf_n(jf)%ln_tint
836         sdf(jf)%ln_clim    = sdf_n(jf)%ln_clim
837         sdf(jf)%cltype     = sdf_n(jf)%cltype
838         sdf(jf)%wgtname = " "
839         IF( LEN( TRIM(sdf_n(jf)%wname) ) > 0 )   sdf(jf)%wgtname = TRIM( cdir )//TRIM( sdf_n(jf)%wname )
840         sdf(jf)%vcomp   = sdf_n(jf)%vcomp
841      END DO
842
843      IF(lwp) THEN      ! control print
844         WRITE(numout,*)
845         WRITE(numout,*) TRIM( cdcaller )//' : '//TRIM( cdtitle )
846         WRITE(numout,*) (/ ('~', jf = 1, LEN_TRIM( cdcaller ) ) /)
847         WRITE(numout,*) '          '//TRIM( cdnam )//' Namelist'
848         WRITE(numout,*) '          list of files and frequency (>0: in hours ; <0 in months)'
849         DO jf = 1, SIZE(sdf)
850            WRITE(numout,*) '               root filename: '  , TRIM( sdf(jf)%clrootname ),   &
851               &                          ' variable name: '  , TRIM( sdf(jf)%clvar      )
852            WRITE(numout,*) '               frequency: '      ,       sdf(jf)%nfreqh      ,   &
853               &                          ' time interp: '    ,       sdf(jf)%ln_tint     ,   &
854               &                          ' climatology: '    ,       sdf(jf)%ln_clim     ,   &
855               &                          ' weights    : '    , TRIM( sdf(jf)%wgtname    ),   &
856               &                          ' pairing    : '    , TRIM( sdf(jf)%vcomp      ),   &
857               &                          ' data type: '      ,       sdf(jf)%cltype
858            call flush(numout)
859         END DO
860      ENDIF
861     
862   END SUBROUTINE fld_fill
863
864
865   SUBROUTINE wgt_list( sd, kwgt )
866      !!---------------------------------------------------------------------
867      !!                    ***  ROUTINE wgt_list  ***
868      !!
869      !! ** Purpose :   search array of WGTs and find a weights file
870      !!                entry, or return a new one adding it to the end
871      !!                if it is a new entry, the weights data is read in and
872      !!                restructured (fld_weight)
873      !!----------------------------------------------------------------------
874      TYPE( FLD ), INTENT(in   ) ::   sd        ! field with name of weights file
875      INTEGER    , INTENT(inout) ::   kwgt      ! index of weights
876      !!
877      INTEGER ::   kw, nestid   ! local integer
878      LOGICAL ::   found        ! local logical
879      !!----------------------------------------------------------------------
880      !
881      !! search down linked list
882      !! weights filename is either present or we hit the end of the list
883      found = .FALSE.
884
885      !! because agrif nest part of filenames are now added in iom_open
886      !! to distinguish between weights files on the different grids, need to track
887      !! nest number explicitly
888      nestid = 0
889#if defined key_agrif
890      nestid = Agrif_Fixed()
891#endif
892      DO kw = 1, nxt_wgt-1
893         IF( TRIM(ref_wgts(kw)%wgtname) == TRIM(sd%wgtname) .AND. &
894             ref_wgts(kw)%nestid == nestid) THEN
895            kwgt = kw
896            found = .TRUE.
897            EXIT
898         ENDIF
899      END DO
900      IF( .NOT.found ) THEN
901         kwgt = nxt_wgt
902         CALL fld_weight( sd )
903      ENDIF
904      !
905   END SUBROUTINE wgt_list
906
907
908   SUBROUTINE wgt_print( )
909      !!---------------------------------------------------------------------
910      !!                    ***  ROUTINE wgt_print  ***
911      !!
912      !! ** Purpose :   print the list of known weights
913      !!----------------------------------------------------------------------
914      INTEGER ::   kw   !
915      !!----------------------------------------------------------------------
916      !
917      DO kw = 1, nxt_wgt-1
918         WRITE(numout,*) 'weight file:  ',TRIM(ref_wgts(kw)%wgtname)
919         WRITE(numout,*) '      ddims:  ',ref_wgts(kw)%ddims(1),ref_wgts(kw)%ddims(2)
920         WRITE(numout,*) '     numwgt:  ',ref_wgts(kw)%numwgt
921         WRITE(numout,*) '     jpiwgt:  ',ref_wgts(kw)%jpiwgt
922         WRITE(numout,*) '     jpjwgt:  ',ref_wgts(kw)%jpjwgt
923         WRITE(numout,*) '    botleft:  ',ref_wgts(kw)%botleft
924         WRITE(numout,*) '   topright:  ',ref_wgts(kw)%topright
925         IF( ref_wgts(kw)%cyclic ) THEN
926            WRITE(numout,*) '       cyclical'
927            IF( ref_wgts(kw)%overlap > 0 ) WRITE(numout,*) '              with overlap of ', ref_wgts(kw)%overlap
928         ELSE
929            WRITE(numout,*) '       not cyclical'
930         ENDIF
931         IF( ASSOCIATED(ref_wgts(kw)%data_wgt) )  WRITE(numout,*) '       allocated'
932      END DO
933      !
934   END SUBROUTINE wgt_print
935
936
937   SUBROUTINE fld_weight( sd )
938      !!---------------------------------------------------------------------
939      !!                    ***  ROUTINE fld_weight  ***
940      !!
941      !! ** Purpose :   create a new WGT structure and fill in data from 
942      !!                file, restructuring as required
943      !!----------------------------------------------------------------------
944      TYPE( FLD ), INTENT(in) ::   sd   ! field with name of weights file
945      !!
946      INTEGER                           ::   jn            ! dummy loop indices
947      INTEGER                           ::   inum          ! temporary logical unit
948      INTEGER                           ::   id            ! temporary variable id
949      INTEGER                           ::   ipk           ! temporary vertical dimension
950      CHARACTER (len=5)                 ::   aname
951      INTEGER , DIMENSION(3)            ::   ddims
952      INTEGER , POINTER, DIMENSION(:,:) ::   data_src
953      REAL(wp), POINTER, DIMENSION(:,:) ::   data_tmp
954      LOGICAL                           ::   cyclical
955      INTEGER                           ::   zwrap      ! local integer
956      !!----------------------------------------------------------------------
957      !
958      CALL wrk_alloc( jpi,jpj, data_src )   ! integer
959      CALL wrk_alloc( jpi,jpj, data_tmp )
960      !
961      IF( nxt_wgt > tot_wgts ) THEN
962        CALL ctl_stop("fld_weight: weights array size exceeded, increase tot_wgts")
963      ENDIF
964      !
965      !! new weights file entry, add in extra information
966      !! a weights file represents a 2D grid of a certain shape, so we assume that the current
967      !! input data file is representative of all other files to be opened and processed with the
968      !! current weights file
969
970      !! open input data file (non-model grid)
971      CALL iom_open( sd%clname, inum, ldiof =  LEN(TRIM(sd%wgtname)) > 0 )
972
973      !! get dimensions
974      id = iom_varid( inum, sd%clvar, ddims )
975
976      !! close it
977      CALL iom_close( inum )
978
979      !! now open the weights file
980
981      CALL iom_open ( sd%wgtname, inum )   ! interpolation weights
982      IF ( inum > 0 ) THEN
983
984         !! determine whether we have an east-west cyclic grid
985         !! from global attribute called "ew_wrap" in the weights file
986         !! note that if not found, iom_getatt returns -999 and cyclic with no overlap is assumed
987         !! since this is the most common forcing configuration
988
989         CALL iom_getatt(inum, 'ew_wrap', zwrap)
990         IF( zwrap >= 0 ) THEN
991            cyclical = .TRUE.
992         ELSE IF( zwrap == -999 ) THEN
993            cyclical = .TRUE.
994            zwrap = 0
995         ELSE
996            cyclical = .FALSE.
997         ENDIF
998
999         ref_wgts(nxt_wgt)%ddims(1) = ddims(1)
1000         ref_wgts(nxt_wgt)%ddims(2) = ddims(2)
1001         ref_wgts(nxt_wgt)%wgtname = sd%wgtname
1002         ref_wgts(nxt_wgt)%overlap = zwrap
1003         ref_wgts(nxt_wgt)%cyclic = cyclical
1004         ref_wgts(nxt_wgt)%nestid = 0
1005#if defined key_agrif
1006         ref_wgts(nxt_wgt)%nestid = Agrif_Fixed()
1007#endif
1008         !! weights file is stored as a set of weights (wgt01->wgt04 or wgt01->wgt16)
1009         !! for each weight wgtNN there is an integer array srcNN which gives the point in
1010         !! the input data grid which is to be multiplied by the weight
1011         !! they are both arrays on the model grid so the result of the multiplication is
1012         !! added into an output array on the model grid as a running sum
1013
1014         !! two possible cases: bilinear (4 weights) or bicubic (16 weights)
1015         id = iom_varid(inum, 'src05', ldstop=.FALSE.)
1016         IF( id <= 0) THEN
1017            ref_wgts(nxt_wgt)%numwgt = 4
1018         ELSE
1019            ref_wgts(nxt_wgt)%numwgt = 16
1020         ENDIF
1021
1022         ALLOCATE( ref_wgts(nxt_wgt)%data_jpi(jpi,jpj,4) )
1023         ALLOCATE( ref_wgts(nxt_wgt)%data_jpj(jpi,jpj,4) )
1024         ALLOCATE( ref_wgts(nxt_wgt)%data_wgt(jpi,jpj,ref_wgts(nxt_wgt)%numwgt) )
1025
1026         DO jn = 1,4
1027            aname = ' '
1028            WRITE(aname,'(a3,i2.2)') 'src',jn
1029            data_tmp(:,:) = 0
1030            CALL iom_get ( inum, jpdom_data, aname, data_tmp(:,:) )
1031            data_src(:,:) = INT(data_tmp(:,:))
1032            ref_wgts(nxt_wgt)%data_jpj(:,:,jn) = 1 + (data_src(:,:)-1) / ref_wgts(nxt_wgt)%ddims(1)
1033            ref_wgts(nxt_wgt)%data_jpi(:,:,jn) = data_src(:,:) - ref_wgts(nxt_wgt)%ddims(1)*(ref_wgts(nxt_wgt)%data_jpj(:,:,jn)-1)
1034         END DO
1035
1036         DO jn = 1, ref_wgts(nxt_wgt)%numwgt
1037            aname = ' '
1038            WRITE(aname,'(a3,i2.2)') 'wgt',jn
1039            ref_wgts(nxt_wgt)%data_wgt(:,:,jn) = 0.0
1040            CALL iom_get ( inum, jpdom_data, aname, ref_wgts(nxt_wgt)%data_wgt(:,:,jn) )
1041         END DO
1042         CALL iom_close (inum)
1043 
1044         ! find min and max indices in grid
1045         ref_wgts(nxt_wgt)%botleft(1) = MINVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1046         ref_wgts(nxt_wgt)%botleft(2) = MINVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1047         ref_wgts(nxt_wgt)%topright(1) = MAXVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1048         ref_wgts(nxt_wgt)%topright(2) = MAXVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1049
1050         ! and therefore dimensions of the input box
1051         ref_wgts(nxt_wgt)%jpiwgt = ref_wgts(nxt_wgt)%topright(1) - ref_wgts(nxt_wgt)%botleft(1) + 1
1052         ref_wgts(nxt_wgt)%jpjwgt = ref_wgts(nxt_wgt)%topright(2) - ref_wgts(nxt_wgt)%botleft(2) + 1
1053
1054         ! shift indexing of source grid
1055         ref_wgts(nxt_wgt)%data_jpi(:,:,:) = ref_wgts(nxt_wgt)%data_jpi(:,:,:) - ref_wgts(nxt_wgt)%botleft(1) + 1
1056         ref_wgts(nxt_wgt)%data_jpj(:,:,:) = ref_wgts(nxt_wgt)%data_jpj(:,:,:) - ref_wgts(nxt_wgt)%botleft(2) + 1
1057
1058         ! create input grid, give it a halo to allow gradient calculations
1059         ! SA: +3 stencil is a patch to avoid out-of-bound computation in some configuration.
1060         ! a more robust solution will be given in next release
1061         ipk =  SIZE(sd%fnow, 3)
1062         ALLOCATE( ref_wgts(nxt_wgt)%fly_dta(ref_wgts(nxt_wgt)%jpiwgt+3, ref_wgts(nxt_wgt)%jpjwgt+3 ,ipk) )
1063         IF( ref_wgts(nxt_wgt)%cyclic ) ALLOCATE( ref_wgts(nxt_wgt)%col(1,ref_wgts(nxt_wgt)%jpjwgt+3,ipk) )
1064
1065         nxt_wgt = nxt_wgt + 1
1066
1067      ELSE
1068         CALL ctl_stop( '    fld_weight : unable to read the file ' )
1069      ENDIF
1070
1071      CALL wrk_dealloc( jpi,jpj, data_src )   ! integer
1072      CALL wrk_dealloc( jpi,jpj, data_tmp )
1073      !
1074   END SUBROUTINE fld_weight
1075
1076
1077   SUBROUTINE fld_interp( num, clvar, kw, kk, dta, nrec )
1078      !!---------------------------------------------------------------------
1079      !!                    ***  ROUTINE fld_interp  ***
1080      !!
1081      !! ** Purpose :   apply weights to input gridded data to create data
1082      !!                on model grid
1083      !!----------------------------------------------------------------------
1084      INTEGER                   , INTENT(in   ) ::   num     ! stream number
1085      CHARACTER(LEN=*)          , INTENT(in   ) ::   clvar   ! variable name
1086      INTEGER                   , INTENT(in   ) ::   kw      ! weights number
1087      INTEGER                   , INTENT(in   ) ::   kk      ! vertical dimension of kk
1088      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   dta     ! output field on model grid
1089      INTEGER                   , INTENT(in   ) ::   nrec    ! record number to read (ie time slice)
1090      !!
1091      INTEGER, DIMENSION(3) ::   rec1,recn   ! temporary arrays for start and length
1092      INTEGER ::  jk, jn, jm           ! loop counters
1093      INTEGER ::  ni, nj               ! lengths
1094      INTEGER ::  jpimin,jpiwid        ! temporary indices
1095      INTEGER ::  jpjmin,jpjwid        ! temporary indices
1096      INTEGER ::  jpi1,jpi2,jpj1,jpj2  ! temporary indices
1097      !!----------------------------------------------------------------------
1098      !
1099      !! for weighted interpolation we have weights at four corners of a box surrounding
1100      !! a model grid point, each weight is multiplied by a grid value (bilinear case)
1101      !! or by a grid value and gradients at the corner point (bicubic case)
1102      !! so we need to have a 4 by 4 subgrid surrounding each model point to cover both cases
1103
1104      !! sub grid from non-model input grid which encloses all grid points in this nemo process
1105      jpimin = ref_wgts(kw)%botleft(1)
1106      jpjmin = ref_wgts(kw)%botleft(2)
1107      jpiwid = ref_wgts(kw)%jpiwgt
1108      jpjwid = ref_wgts(kw)%jpjwgt
1109
1110      !! when reading in, expand this sub-grid by one halo point all the way round for calculating gradients
1111      rec1(1) = MAX( jpimin-1, 1 )
1112      rec1(2) = MAX( jpjmin-1, 1 )
1113      rec1(3) = 1
1114      recn(1) = MIN( jpiwid+2, ref_wgts(kw)%ddims(1)-rec1(1)+1 )
1115      recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
1116      recn(3) = kk
1117
1118      !! where we need to put it in the non-nemo grid fly_dta
1119      !! note that jpi1 and jpj1 only differ from 1 when jpimin and jpjmin are 1
1120      !! (ie at the extreme west or south of the whole input grid) and similarly for jpi2 and jpj2
1121      jpi1 = 2 + rec1(1) - jpimin
1122      jpj1 = 2 + rec1(2) - jpjmin
1123      jpi2 = jpi1 + recn(1) - 1
1124      jpj2 = jpj1 + recn(2) - 1
1125
1126      ref_wgts(kw)%fly_dta(:,:,:) = 0.0
1127      SELECT CASE( SIZE(ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:),3) )
1128      CASE(1)
1129           CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,1), nrec, rec1, recn)
1130      CASE DEFAULT
1131           CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:), nrec, rec1, recn)
1132      END SELECT 
1133
1134      !! first four weights common to both bilinear and bicubic
1135      !! data_jpi, data_jpj have already been shifted to (1,1) corresponding to botleft
1136      !! note that we have to offset by 1 into fly_dta array because of halo
1137      dta(:,:,:) = 0.0
1138      DO jk = 1,4
1139        DO jn = 1, jpj
1140          DO jm = 1,jpi
1141            ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1142            nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1143            dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk) * ref_wgts(kw)%fly_dta(ni+1,nj+1,:)
1144          END DO
1145        END DO
1146      END DO
1147
1148      IF (ref_wgts(kw)%numwgt .EQ. 16) THEN
1149
1150        !! fix up halo points that we couldnt read from file
1151        IF( jpi1 == 2 ) THEN
1152           ref_wgts(kw)%fly_dta(jpi1-1,:,:) = ref_wgts(kw)%fly_dta(jpi1,:,:)
1153        ENDIF
1154        IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
1155           ref_wgts(kw)%fly_dta(jpi2+1,:,:) = ref_wgts(kw)%fly_dta(jpi2,:,:)
1156        ENDIF
1157        IF( jpj1 == 2 ) THEN
1158           ref_wgts(kw)%fly_dta(:,jpj1-1,:) = ref_wgts(kw)%fly_dta(:,jpj1,:)
1159        ENDIF
1160        IF( jpj2 + jpjmin - 1 == ref_wgts(kw)%ddims(2)+1 .AND. jpj2 .lt. jpjwid+2 ) THEN
1161           ref_wgts(kw)%fly_dta(:,jpj2+1,:) = 2.0*ref_wgts(kw)%fly_dta(:,jpj2,:) - ref_wgts(kw)%fly_dta(:,jpj2-1,:)
1162        ENDIF
1163
1164        !! if data grid is cyclic we can do better on east-west edges
1165        !! but have to allow for whether first and last columns are coincident
1166        IF( ref_wgts(kw)%cyclic ) THEN
1167           rec1(2) = MAX( jpjmin-1, 1 )
1168           recn(1) = 1
1169           recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
1170           jpj1 = 2 + rec1(2) - jpjmin
1171           jpj2 = jpj1 + recn(2) - 1
1172           IF( jpi1 == 2 ) THEN
1173              rec1(1) = ref_wgts(kw)%ddims(1) - ref_wgts(kw)%overlap
1174              SELECT CASE( SIZE( ref_wgts(kw)%col(:,jpj1:jpj2,:),3) )
1175              CASE(1)
1176                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,1), nrec, rec1, recn)
1177              CASE DEFAULT
1178                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, rec1, recn)
1179              END SELECT     
1180              ref_wgts(kw)%fly_dta(jpi1-1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
1181           ENDIF
1182           IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
1183              rec1(1) = 1 + ref_wgts(kw)%overlap
1184              SELECT CASE( SIZE( ref_wgts(kw)%col(:,jpj1:jpj2,:),3) )
1185              CASE(1)
1186                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,1), nrec, rec1, recn)
1187              CASE DEFAULT
1188                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, rec1, recn)
1189              END SELECT
1190              ref_wgts(kw)%fly_dta(jpi2+1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
1191           ENDIF
1192        ENDIF
1193
1194        ! gradient in the i direction
1195        DO jk = 1,4
1196          DO jn = 1, jpj
1197            DO jm = 1,jpi
1198              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1199              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1200              dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+4) * 0.5 *         &
1201                               (ref_wgts(kw)%fly_dta(ni+2,nj+1,:) - ref_wgts(kw)%fly_dta(ni,nj+1,:))
1202            END DO
1203          END DO
1204        END DO
1205
1206        ! gradient in the j direction
1207        DO jk = 1,4
1208          DO jn = 1, jpj
1209            DO jm = 1,jpi
1210              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1211              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1212              dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+8) * 0.5 *         &
1213                               (ref_wgts(kw)%fly_dta(ni+1,nj+2,:) - ref_wgts(kw)%fly_dta(ni+1,nj,:))
1214            END DO
1215          END DO
1216        END DO
1217
1218         ! gradient in the ij direction
1219         DO jk = 1,4
1220            DO jn = 1, jpj
1221               DO jm = 1,jpi
1222                  ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1223                  nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1224                  dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+12) * 0.25 * ( &
1225                               (ref_wgts(kw)%fly_dta(ni+2,nj+2,:) - ref_wgts(kw)%fly_dta(ni  ,nj+2,:)) -   &
1226                               (ref_wgts(kw)%fly_dta(ni+2,nj  ,:) - ref_wgts(kw)%fly_dta(ni  ,nj  ,:)))
1227               END DO
1228            END DO
1229         END DO
1230         !
1231      END IF
1232      !
1233   END SUBROUTINE fld_interp
1234
1235
1236   FUNCTION ksec_week( cdday )
1237      !!---------------------------------------------------------------------
1238      !!                    ***  FUNCTION kshift_week ***
1239      !!
1240      !! ** Purpose : 
1241      !!---------------------------------------------------------------------
1242      CHARACTER(len=*), INTENT(in)   ::   cdday   !3 first letters of the first day of the weekly file
1243      !!
1244      INTEGER                        ::   ksec_week  ! output variable
1245      INTEGER                        ::   ijul       !temp variable
1246      INTEGER                        ::   ishift     !temp variable
1247      CHARACTER(len=3),DIMENSION(7)  ::   cl_week 
1248      !!----------------------------------------------------------------------
1249      cl_week = (/"sun","sat","fri","thu","wed","tue","mon"/)
1250      DO ijul = 1, 7
1251         IF( cl_week(ijul) == TRIM(cdday) ) EXIT
1252      END DO
1253      IF( ijul .GT. 7 )   CALL ctl_stop( 'ksec_week: wrong day for sdjf%cltype(6:8): '//TRIM(cdday) )
1254      !
1255      ishift = ijul * NINT(rday)
1256      !
1257      ksec_week = nsec_week + ishift
1258      ksec_week = MOD( ksec_week, 7*NINT(rday) )
1259      !
1260   END FUNCTION ksec_week
1261
1262   !!======================================================================
1263END MODULE fldread
Note: See TracBrowser for help on using the repository browser.