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/NERC/dev_r3874_FASTNEt/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/NERC/dev_r3874_FASTNEt/NEMOGCM/NEMO/OPA_SRC/SBC/fldread.F90 @ 6736

Last change on this file since 6736 was 6736, checked in by jamesharle, 8 years ago

FASTNEt code modifications

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