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

source: trunk/NEMOGCM/NEMO/OPA_SRC/SBC/fldread.F90 @ 3294

Last change on this file since 3294 was 3294, checked in by rblod, 12 years ago

Merge of 3.4beta into the trunk

  • Property svn:keywords set to Id
File size: 67.4 KB
RevLine 
[888]1MODULE fldread
2   !!======================================================================
3   !!                       ***  MODULE  fldread  ***
4   !! Ocean forcing:  read input field for surface boundary condition
5   !!=====================================================================
[2715]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
[888]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
[1275]20   USE geo2ocean       ! for vector rotation on to model grid
[2715]21   USE lib_mpp         ! MPP library
[3294]22   USE wrk_nemo        ! work arrays
[2715]23   USE ioipsl, ONLY :   ymds2ju, ju2ymds   ! for calendar
[888]24
25   IMPLICIT NONE
26   PRIVATE   
[3294]27 
28   PUBLIC   fld_map    ! routine called by tides_init
[888]29
30   TYPE, PUBLIC ::   FLD_N      !: Namelist field informations
[1730]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)
[2528]36      CHARACTER(len = 8)   ::   cltype      ! type of data file 'daily', 'monthly' or yearly'
[1730]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
[2715]39      !                                     ! a string starting with "U" or "V" for each component   
40      !                                     ! chars 2 onwards identify which components go together 
[888]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
[1730]46      INTEGER                         ::   nfreqh       ! frequency of each flux file
[888]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)
[1132]49      LOGICAL                         ::   ln_clim      ! climatology or not (T/F)
[2528]50      CHARACTER(len = 8)              ::   cltype       ! type of data file 'daily', 'monthly' or yearly'
[1132]51      INTEGER                         ::   num          ! iom id of the jpfld files to be read
[1730]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)
[2528]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
[1275]56      CHARACTER(len = 256)            ::   wgtname      ! current name of the NetCDF weight file acting as a key
[2715]57      !                                                 ! into the WGTLIST structure
[1275]58      CHARACTER(len = 34)             ::   vcomp        ! symbolic name for a vector component that needs rotation
[2528]59      LOGICAL                         ::   rotn         ! flag to indicate whether field has been rotated
[888]60   END TYPE FLD
61
[3294]62   TYPE, PUBLIC ::   MAP_POINTER      !: Array of integer pointers to 1D arrays
63      INTEGER, POINTER   ::  ptr(:)
64   END TYPE MAP_POINTER
65
[1275]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
[2715]75      !                                                         ! current processor grid
[1275]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
[2528]81      INTEGER                                 ::   overlap      ! =0 when cyclic grid has no overlapping EW columns
[2715]82      !                                                         ! =>1 when they have one or more overlapping columns     
83      !                                                         ! =-1 not cyclic
[1275]84      LOGICAL                                 ::   cyclic       ! east-west cyclic or not
[2528]85      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpi     ! array of source integers
86      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpj     ! array of source integers
[1275]87      REAL(wp), DIMENSION(:,:,:), POINTER     ::   data_wgt     ! array of weights on model grid
[2528]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
[1275]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
[1132]98   PUBLIC   fld_read, fld_fill   ! called by sbc... modules
[888]99
100   !!----------------------------------------------------------------------
[2528]101   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[1156]102   !! $Id$
[2715]103   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[888]104   !!----------------------------------------------------------------------
105CONTAINS
106
[3294]107   SUBROUTINE fld_read( kt, kn_fsbc, sd, map, jit, time_offset )
[888]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
[1132]120      INTEGER  , INTENT(in   )               ::   kn_fsbc   ! sbc computation period (in time step)
[888]121      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
[3294]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.
[888]128      !!
[2528]129      INTEGER  ::   imf        ! size of the structure sd
[1132]130      INTEGER  ::   jf         ! dummy indices
[1730]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
[2323]133      INTEGER  ::   isecsbc    ! number of seconds between Jan. 1st 00h of nit000 year and the middle of sbc time step
[3294]134      INTEGER  ::   itime_add  ! local time offset variable
[1628]135      LOGICAL  ::   llnxtyr    ! open next year  file?
136      LOGICAL  ::   llnxtmth   ! open next month file?
137      LOGICAL  ::   llstop     ! stop is the file does not exist
[3294]138      LOGICAL  ::   ll_firstcall ! true if this is the first call to fld_read for this set of fields
[1132]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
[1191]141      CHARACTER(LEN=1000) ::   clfmt   ! write format
[888]142      !!---------------------------------------------------------------------
[3294]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         
[2528]153      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
[3294]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
[1275]160      imf = SIZE( sd )
[2323]161      !
[3294]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
[2528]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         !                                         ! ====================================== !
[888]178         !
[2528]179         DO jf = 1, imf                            ! ---   loop over field   --- !
180           
[3294]181            IF( isecsbc > sd(jf)%nrec_a(2) .OR. ll_firstcall ) THEN  ! read/update the after data?
[888]182
[2528]183               IF( sd(jf)%ln_tint ) THEN                             ! swap before record field and informations
184                  sd(jf)%nrec_b(:) = sd(jf)%nrec_a(:)
[888]185!CDIR COLLAPSE
[2528]186                  sd(jf)%fdta(:,:,:,1) = sd(jf)%fdta(:,:,:,2)
187               ENDIF
[1132]188
[3294]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
[1132]194
[2528]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
[1132]200
[2528]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
[1132]213                  ENDIF
214
[2528]215                  ! do we need next file data?
216                  IF( sd(jf)%nrec_a(1) > ireclast ) THEN
[1132]217
[2528]218                     sd(jf)%nrec_a(1) = 1              ! force to read the first record of the next file
[1628]219
[2528]220                     IF( .NOT. sd(jf)%ln_clim ) THEN   ! close the current file and open a new one.
[1628]221
[2528]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?
[1132]224
[2528]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
[1132]230
[2528]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
[1132]242                     ENDIF
[2528]243                  ENDIF
[1132]244
[2528]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 )
[1132]251               ENDIF
252
[2528]253               ! read after data
[3294]254               IF( PRESENT(map) ) THEN
255                  CALL fld_get( sd(jf), map(jf)%ptr )
256               ELSE
257                  CALL fld_get( sd(jf) )
258               ENDIF
[2528]259
[1275]260            ENDIF
[2528]261         END DO                                    ! --- end loop over field --- !
[1132]262
[2528]263         CALL fld_rot( kt, sd )                    ! rotate vector fiels if needed
[888]264
[2528]265         DO jf = 1, imf                            ! ---   loop over field   --- !
[888]266            !
[2528]267            IF( sd(jf)%ln_tint ) THEN              ! temporal interpolation
[1191]268               IF(lwp .AND. kt - nit000 <= 100 ) THEN
[2528]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, ')')"
[3294]271                  WRITE(numout, clfmt)  TRIM( sd(jf)%clvar ), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,   &           
[1730]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
[3294]273                  WRITE(numout, *) 'itime_add is : ',itime_add
[1191]274               ENDIF
[2528]275               ! temporal interpolation weights
[2323]276               ztinta =  REAL( isecsbc - sd(jf)%nrec_b(2), wp ) / REAL( sd(jf)%nrec_a(2) - sd(jf)%nrec_b(2), wp )
[1132]277               ztintb =  1. - ztinta
[888]278!CDIR COLLAPSE
[2528]279               sd(jf)%fnow(:,:,:) = ztintb * sd(jf)%fdta(:,:,:,1) + ztinta * sd(jf)%fdta(:,:,:,2)
280            ELSE   ! nothing to do...
[1191]281               IF(lwp .AND. kt - nit000 <= 100 ) THEN
[2528]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
[1191]286               ENDIF
[888]287            ENDIF
288            !
[2528]289            IF( kt == nitend - kn_fsbc + 1 )   CALL iom_close( sd(jf)%num )   ! Close the input files
[1132]290
[2528]291         END DO                                    ! --- end loop over field --- !
292         !
293         !                                         ! ====================================== !
294      ENDIF                                        ! update field at each kn_fsbc time-step !
295      !                                            ! ====================================== !
296      !
[888]297   END SUBROUTINE fld_read
298
299
[3294]300   SUBROUTINE fld_init( kn_fsbc, sdjf, map )
[888]301      !!---------------------------------------------------------------------
[1132]302      !!                    ***  ROUTINE fld_init  ***
303      !!
304      !! ** Purpose :  - if time interpolation, read before data
305      !!               - open current year file
306      !!----------------------------------------------------------------------
[2528]307      INTEGER  , INTENT(in   ) ::   kn_fsbc   ! sbc computation period (in time step)
308      TYPE(FLD), INTENT(inout) ::   sdjf      ! input field related variables
[3294]309      INTEGER  , INTENT(in), OPTIONAL, DIMENSION(:) :: map ! global-to-local mapping indices
[1132]310      !!
[2528]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
[1191]320      CHARACTER(LEN=1000) ::   clfmt   ! write format
[1132]321      !!---------------------------------------------------------------------
[2528]322     
[1132]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
[2528]326      llprevyr   = .FALSE.
327      llprevmth  = .FALSE.
328      llprevweek = .FALSE.
329      llprevday  = .FALSE.
330      isec_week  = 0
[1132]331           
[2528]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
[1132]337      ! define record informations
[2528]338      CALL fld_rec( kn_fsbc, sdjf, ldbefore = .TRUE. )  ! return before values in sdjf%nrec_a (as we will swap it later)
[1132]339
[2528]340      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
[2323]341
[1132]342      IF( sdjf%ln_tint ) THEN ! we need to read the previous record and we will put it in the current record structure
[2528]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?
[1628]356                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
[2528]357               ELSE                                           ! yearly file
358                  sdjf%nrec_a(1) = 12                                                      ! force to read december mean
[1628]359                  llprevyr = .NOT. sdjf%ln_clim                                            ! use previous year  file?
360               ENDIF
[2528]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?
[1628]365                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
[2528]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?
[1628]373                  llprevmth = llprevday .AND. nday   == 1                                  ! use previous month file?
374                  llprevyr  = llprevmth .AND. nmonth == 1                                  ! use previous year  file?
[2528]375               ELSE                                           ! yearly file
376                  sdjf%nrec_a(1) = 24 * nyear_len(0) / sdjf%nfreqh                         ! last record of previous year
[1628]377                  llprevyr = .NOT. sdjf%ln_clim                                            ! use previous year  file?
[1132]378               ENDIF
379            ENDIF
380         ENDIF
[2528]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 )
[1132]393
[1628]394         ! if previous year/month/day file does not exist, we switch to the current year/month/day
[1818]395         IF( llprev .AND. sdjf%num <= 0 ) THEN
[2528]396            CALL ctl_warn( 'previous year/month/week/day file: '//TRIM(sdjf%clname)//   &
397               &           ' not present -> back to current year/month/week/day' )
[1628]398            ! we force to read the first record of the current year/month/day instead of last record of previous year/month/day
[2528]399            llprev = .FALSE.
400            sdjf%nrec_a(1) = 1
[1628]401            CALL fld_clopn( sdjf, nyear, nmonth, nday )
[1132]402         ENDIF
403         
[1730]404         IF( llprev ) THEN   ! check if the last record sdjf%nrec_n(1) exists in the file
[1132]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
[2528]408            sdjf%nrec_a(1) = MIN( sdjf%nrec_a(1), inrec )   ! make sure we select an existing record
[1132]409         ENDIF
410
[2528]411         ! read before data
[3294]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
[1132]417
[1191]418         clfmt = "('fld_init : time-interpolation for ', a, ' read previous record = ', i4, ' at time = ', f7.2, ' days')"
[2528]419         IF(lwp) WRITE(numout, clfmt) TRIM(sdjf%clvar), sdjf%nrec_a(1), REAL(sdjf%nrec_a(2),wp)/rday
[1132]420
[2528]421         IF( llprev )   CALL iom_close( sdjf%num )          ! force to close previous year file (-> redefine sdjf%num to 0)
[1132]422
423      ENDIF
424
[2528]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 
[2715]444      !
[1132]445   END SUBROUTINE fld_init
446
447
[3294]448   SUBROUTINE fld_rec( kn_fsbc, sdjf, ldbefore, jit, time_offset )
[1132]449      !!---------------------------------------------------------------------
[888]450      !!                    ***  ROUTINE fld_rec  ***
451      !!
[2528]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)
[888]458      !!----------------------------------------------------------------------
[2528]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.)
[3294]462      INTEGER  , INTENT(in   ), OPTIONAL ::   jit       ! index of barotropic subcycle
[2528]463                                                        ! used only if sdjf%ln_tint = .TRUE.
[3294]464      INTEGER  , INTENT(in   ), OPTIONAL ::   time_offset  ! Offset of required time level compared to "now"
465                                                           ! time level in units of time steps.
[888]466      !!
[2528]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
[3294]472      INTEGER  ::   itime_add   ! local time offset variable
[1132]473      REAL(wp) ::   ztmp        ! temporary variable
[888]474      !!----------------------------------------------------------------------
475      !
[2528]476      ! Note that shifting time to be centrered in the middle of sbc time step impacts only nsec_* variables of the calendar
[2323]477      !
[2528]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      !
[3294]482      itime_add = 0
483      IF( PRESENT(time_offset) ) itime_add = time_offset
484      !
[2528]485      !                                      ! =========== !
486      IF    ( sdjf%nfreqh == -12 ) THEN      ! yearly mean
487         !                                   ! =========== !
[888]488         !
[1132]489         IF( sdjf%ln_tint ) THEN                 ! time interpolation, shift by 1/2 record
490            !
491            !                  INT( ztmp )
492            !                     /|\
493            !                    1 |    *----
494            !                    0 |----(             
495            !                      |----+----|--> time
[2528]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
[3294]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
[2528]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
[1132]529            !                      0   /|\   1   (nday/nmonth_len(nmonth))
530            !                           |   
531            !                           |   
532            !       forcing record :  nmonth
533            !                           
[2528]534            ztmp = REAL( nday, wp ) / REAL( nmonth_len(nmonth), wp ) + 0.5
[3294]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
[2528]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)
[888]551         ENDIF
552         !
[2528]553         !                                   ! ================================ !
554      ELSE                                   ! higher frequency mean (in hours)
555         !                                   ! ================================ !
[888]556         !
[2528]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
[1132]559         ! number of second since the beginning of the file
[2528]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
[1132]564         ENDIF
[2323]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
[3294]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
[1132]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
[1730]580            !                      0 /|\ 1 /|\ 2 /|\ 3 (nsec_year/ifreq_sec) or (nsec_month/ifreq_sec)
[1132]581            !                         |     |     |
582            !                         |     |     |
583            !       forcing record :  1     2     3
584            !                   
[2528]585            ztmp= ztmp / REAL(ifreq_sec, wp) + 0.5
586         ELSE                                   ! no time interpolation
[1132]587            !
588            !                  INT( ztmp )
589            !                     /|\
590            !                    2 |           *-----(
591            !                    1 |     *-----(
592            !                    0 |-----(             
593            !                      |--+--|--+--|--+--|--> time
[1730]594            !                      0 /|\ 1 /|\ 2 /|\ 3 (nsec_year/ifreq_sec) or (nsec_month/ifreq_sec)
[1132]595            !                         |     |     |
596            !                         |     |     |
597            !       forcing record :  1     2     3
598            !                           
[2528]599            ztmp= ztmp / REAL(ifreq_sec, wp)
[1132]600         ENDIF
[2528]601         sdjf%nrec_a(1) = 1 + INT( ztmp ) - COUNT((/llbefore/))   ! record nomber to be read
[1132]602
[2528]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
[888]614         !
615      ENDIF
616      !
[1132]617   END SUBROUTINE fld_rec
618
619
[3294]620   SUBROUTINE fld_get( sdjf, map )
[2528]621      !!---------------------------------------------------------------------
[3294]622      !!                    ***  ROUTINE fld_get  ***
[2528]623      !!
624      !! ** Purpose :   read the data
625      !!----------------------------------------------------------------------
[2715]626      TYPE(FLD), INTENT(inout) ::   sdjf   ! input field related variables
[3294]627      INTEGER  , INTENT(in), OPTIONAL, DIMENSION(:) :: map ! global-to-local mapping indices
[2528]628      !!
[2715]629      INTEGER                  ::   ipk    ! number of vertical levels of sdjf%fdta ( 2D: ipk=1 ; 3D: ipk=jpk )
630      INTEGER                  ::   iw     ! index into wgts array
[2528]631      !!---------------------------------------------------------------------
632           
633      ipk = SIZE( sdjf%fnow, 3 )
[3294]634
635      IF( PRESENT(map) ) THEN
636         IF( sdjf%ln_tint ) THEN   ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fdta(:,:,:,2), sdjf%nrec_a(1), map )
637         ELSE                      ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fnow(:,:,:  ), sdjf%nrec_a(1), map )
638         ENDIF
639      ELSE IF( LEN(TRIM(sdjf%wgtname)) > 0 ) THEN
[2528]640         CALL wgt_list( sdjf, iw )
641         IF( sdjf%ln_tint ) THEN   ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw , ipk  , sdjf%fdta(:,:,:,2), sdjf%nrec_a(1) )
642         ELSE                      ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw , ipk  , sdjf%fnow(:,:,:  ), sdjf%nrec_a(1) )
643         ENDIF
644      ELSE
645         SELECT CASE( ipk )
646         CASE(1)   
647            IF( sdjf%ln_tint ) THEN   ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fdta(:,:,1,2), sdjf%nrec_a(1) )
648            ELSE                      ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fnow(:,:,1  ), sdjf%nrec_a(1) )
649            ENDIF
650         CASE DEFAULT
651            IF( sdjf%ln_tint ) THEN   ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fdta(:,:,:,2), sdjf%nrec_a(1) )
652            ELSE                      ;   CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fnow(:,:,:  ), sdjf%nrec_a(1) )
653            ENDIF
654         END SELECT
655      ENDIF
656      !
657      sdjf%rotn = .false.   ! vector not yet rotated
658
659   END SUBROUTINE fld_get
660
[3294]661   SUBROUTINE fld_map( num, clvar, dta, nrec, map )
662      !!---------------------------------------------------------------------
663      !!                    ***  ROUTINE fld_get  ***
664      !!
665      !! ** Purpose :   read global data from file and map onto local data
666      !!                using a general mapping (for open boundaries)
667      !!----------------------------------------------------------------------
668#if defined key_bdy
669      USE bdy_oce, ONLY:  dta_global         ! workspace to read in global data arrays
670#endif
[2528]671
[3294]672      INTEGER                   , INTENT(in ) ::   num     ! stream number
673      CHARACTER(LEN=*)          , INTENT(in ) ::   clvar   ! variable name
674      REAL(wp), DIMENSION(:,:,:), INTENT(out) ::   dta   ! output field on model grid (2 dimensional)
675      INTEGER                   , INTENT(in ) ::   nrec    ! record number to read (ie time slice)
676      INTEGER,  DIMENSION(:)    , INTENT(in ) ::   map     ! global-to-local mapping indices
677      !!
678      INTEGER                                 ::   ipi      ! length of boundary data on local process
679      INTEGER                                 ::   ipj      ! length of dummy dimension ( = 1 )
680      INTEGER                                 ::   ipk      ! number of vertical levels of dta ( 2D: ipk=1 ; 3D: ipk=jpk )
681      INTEGER                                 ::   ilendta  ! length of data in file
682      INTEGER                                 ::   idvar    ! variable ID
683      INTEGER                                 ::   ib, ik   ! loop counters
684      INTEGER                                 ::   ierr
685      REAL(wp), POINTER, DIMENSION(:,:,:)     ::   dta_read ! work space for global data
686      !!---------------------------------------------------------------------
687           
688#if defined key_bdy
689      dta_read => dta_global
690#endif
691
692      ipi = SIZE( dta, 1 )
693      ipj = 1
694      ipk = SIZE( dta, 3 )
695
696      idvar   = iom_varid( num, clvar )
697      ilendta = iom_file(num)%dimsz(1,idvar)
698      IF(lwp) WRITE(numout,*) 'Dim size for ',TRIM(clvar),' is ', ilendta
699      IF(lwp) WRITE(numout,*) 'Number of levels for ',TRIM(clvar),' is ', ipk
700
701      SELECT CASE( ipk )
702      CASE(1)   
703         CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1    ), nrec )
704      CASE DEFAULT
705         CALL iom_get ( num, jpdom_unknown, clvar, dta_read(1:ilendta,1:ipj,1:ipk), nrec )
706      END SELECT
707      !
708      DO ib = 1, ipi
709         DO ik = 1, ipk
710            dta(ib,1,ik) =  dta_read(map(ib),1,ik)
711         END DO
712      END DO
713
714   END SUBROUTINE fld_map
715
716
[2528]717   SUBROUTINE fld_rot( kt, sd )
718      !!---------------------------------------------------------------------
[3294]719      !!                    ***  ROUTINE fld_rot  ***
[2528]720      !!
721      !! ** Purpose :   Vector fields may need to be rotated onto the local grid direction
[2715]722      !!----------------------------------------------------------------------
[2528]723      INTEGER  , INTENT(in   )               ::   kt        ! ocean time step
724      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
725      !!
[3294]726      INTEGER                           ::   ju, jv, jk   ! loop indices
727      INTEGER                           ::   imf          ! size of the structure sd
728      INTEGER                           ::   ill          ! character length
729      INTEGER                           ::   iv           ! indice of V component
730      REAL(wp), POINTER, DIMENSION(:,:) ::   utmp, vtmp   ! temporary arrays for vector rotation
731      CHARACTER (LEN=100)               ::   clcomp       ! dummy weight name
[2528]732      !!---------------------------------------------------------------------
[2715]733
[3294]734      CALL wrk_alloc( jpi,jpj, utmp, vtmp )
[2715]735
[2528]736      !! (sga: following code should be modified so that pairs arent searched for each time
737      !
738      imf = SIZE( sd )
739      DO ju = 1, imf
740         ill = LEN_TRIM( sd(ju)%vcomp )
741         IF( ill > 0 .AND. .NOT. sd(ju)%rotn ) THEN   ! find vector rotations required             
742             IF( sd(ju)%vcomp(1:1) == 'U' ) THEN      ! east-west component has symbolic name starting with 'U'
743                ! look for the north-south component which has same symbolic name but with 'U' replaced with 'V'
744                clcomp = 'V' // sd(ju)%vcomp(2:ill)   ! works even if ill == 1
745                iv = -1
746                DO jv = 1, imf
747                  IF( TRIM(sd(jv)%vcomp) == TRIM(clcomp) )   iv = jv
748                END DO
749                IF( iv > 0 ) THEN   ! fields ju and iv are two components which need to be rotated together
750                   DO jk = 1, SIZE( sd(ju)%fnow, 3 )
751                      IF( sd(ju)%ln_tint )THEN
752                         CALL rot_rep( sd(ju)%fdta(:,:,jk,2), sd(iv)%fdta(:,:,jk,2), 'T', 'en->i', utmp(:,:) )
753                         CALL rot_rep( sd(ju)%fdta(:,:,jk,2), sd(iv)%fdta(:,:,jk,2), 'T', 'en->j', vtmp(:,:) )
754                         sd(ju)%fdta(:,:,jk,2) = utmp(:,:)   ;   sd(iv)%fdta(:,:,jk,2) = vtmp(:,:)
755                      ELSE
756                         CALL rot_rep( sd(ju)%fnow(:,:,jk  ), sd(iv)%fnow(:,:,jk  ), 'T', 'en->i', utmp(:,:) )
757                         CALL rot_rep( sd(ju)%fnow(:,:,jk  ), sd(iv)%fnow(:,:,jk  ), 'T', 'en->j', vtmp(:,:) )
758                         sd(ju)%fnow(:,:,jk  ) = utmp(:,:)   ;   sd(iv)%fnow(:,:,jk  ) = vtmp(:,:)
759                      ENDIF
760                   END DO
761                   sd(ju)%rotn = .TRUE.               ! vector was rotated
762                   IF( lwp .AND. kt == nit000 )   WRITE(numout,*)   &
763                      &   'fld_read: vector pair ('//TRIM(sd(ju)%clvar)//', '//TRIM(sd(iv)%clvar)//') rotated on to model grid'
764                ENDIF
765             ENDIF
766          ENDIF
767       END DO
[2715]768      !
[3294]769      CALL wrk_dealloc( jpi,jpj, utmp, vtmp )
[2715]770      !
[2528]771   END SUBROUTINE fld_rot
772
773
[1628]774   SUBROUTINE fld_clopn( sdjf, kyear, kmonth, kday, ldstop )
[1132]775      !!---------------------------------------------------------------------
776      !!                    ***  ROUTINE fld_clopn  ***
777      !!
778      !! ** Purpose :   update the file name and open the file
779      !!----------------------------------------------------------------------
[2715]780      TYPE(FLD)        , INTENT(inout) ::   sdjf     ! input field related variables
781      INTEGER          , INTENT(in   ) ::   kyear    ! year value
782      INTEGER          , INTENT(in   ) ::   kmonth   ! month value
783      INTEGER          , INTENT(in   ) ::   kday     ! day value
784      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldstop   ! stop if open to read a non-existing file (default = .TRUE.)
785      !!----------------------------------------------------------------------
[1132]786
787      IF( sdjf%num /= 0 )   CALL iom_close( sdjf%num )   ! close file if already open
788      ! build the new filename if not climatological data
[2528]789      sdjf%clname=TRIM(sdjf%clrootname)
790      !
791      ! note that sdjf%ln_clim is is only acting on presence of the year in the file
792      IF( .NOT. sdjf%ln_clim ) THEN   
793                                         WRITE(sdjf%clname, '(a,"_y",i4.4)' ) TRIM( sdjf%clrootname ), kyear    ! add year
794         IF( sdjf%cltype /= 'yearly' )   WRITE(sdjf%clname, '(a,"m" ,i2.2)' ) TRIM( sdjf%clname     ), kmonth   ! add month
795      ELSE
796         ! build the new filename if climatological data
797         IF( sdjf%cltype /= 'yearly' )   WRITE(sdjf%clname, '(a,"_m",i2.2)' ) TRIM( sdjf%clrootname ), kmonth   ! add month
[888]798      ENDIF
[2528]799      IF( sdjf%cltype == 'daily' .OR. sdjf%cltype(1:4) == 'week' ) &
800            &                            WRITE(sdjf%clname, '(a,"d" ,i2.2)' ) TRIM( sdjf%clname     ), kday     ! add day
801      !
[1319]802      CALL iom_open( sdjf%clname, sdjf%num, ldstop = ldstop, ldiof =  LEN(TRIM(sdjf%wgtname)) > 0 )
[3294]803     !
[1132]804   END SUBROUTINE fld_clopn
805
806
807   SUBROUTINE fld_fill( sdf, sdf_n, cdir, cdcaller, cdtitle, cdnam )
808      !!---------------------------------------------------------------------
809      !!                    ***  ROUTINE fld_fill  ***
810      !!
811      !! ** Purpose :   fill sdf with sdf_n and control print
812      !!----------------------------------------------------------------------
813      TYPE(FLD)  , DIMENSION(:), INTENT(inout) ::   sdf        ! structure of input fields (file informations, fields read)
814      TYPE(FLD_N), DIMENSION(:), INTENT(in   ) ::   sdf_n      ! array of namelist information structures
815      CHARACTER(len=*)         , INTENT(in   ) ::   cdir       ! Root directory for location of flx files
816      CHARACTER(len=*)         , INTENT(in   ) ::   cdcaller   !
817      CHARACTER(len=*)         , INTENT(in   ) ::   cdtitle    !
818      CHARACTER(len=*)         , INTENT(in   ) ::   cdnam      !
[888]819      !
[1132]820      INTEGER  ::   jf       ! dummy indices
821      !!---------------------------------------------------------------------
[888]822
[1132]823      DO jf = 1, SIZE(sdf)
824         sdf(jf)%clrootname = TRIM( cdir )//TRIM( sdf_n(jf)%clname )
[1730]825         sdf(jf)%nfreqh     = sdf_n(jf)%nfreqh
[1132]826         sdf(jf)%clvar      = sdf_n(jf)%clvar
827         sdf(jf)%ln_tint    = sdf_n(jf)%ln_tint
828         sdf(jf)%ln_clim    = sdf_n(jf)%ln_clim
[2528]829         sdf(jf)%cltype     = sdf_n(jf)%cltype
[1275]830         sdf(jf)%wgtname = " "
[1730]831         IF( LEN( TRIM(sdf_n(jf)%wname) ) > 0 )   sdf(jf)%wgtname = TRIM( cdir )//TRIM( sdf_n(jf)%wname )
[1275]832         sdf(jf)%vcomp   = sdf_n(jf)%vcomp
[1132]833      END DO
834
835      IF(lwp) THEN      ! control print
836         WRITE(numout,*)
837         WRITE(numout,*) TRIM( cdcaller )//' : '//TRIM( cdtitle )
838         WRITE(numout,*) (/ ('~', jf = 1, LEN_TRIM( cdcaller ) ) /)
839         WRITE(numout,*) '          '//TRIM( cdnam )//' Namelist'
840         WRITE(numout,*) '          list of files and frequency (>0: in hours ; <0 in months)'
841         DO jf = 1, SIZE(sdf)
842            WRITE(numout,*) '               root filename: '  , TRIM( sdf(jf)%clrootname ),   &
843               &                          ' variable name: '  , TRIM( sdf(jf)%clvar      )
[1730]844            WRITE(numout,*) '               frequency: '      ,       sdf(jf)%nfreqh      ,   &
[1132]845               &                          ' time interp: '    ,       sdf(jf)%ln_tint     ,   &
846               &                          ' climatology: '    ,       sdf(jf)%ln_clim     ,   &
[1275]847               &                          ' weights    : '    , TRIM( sdf(jf)%wgtname    ),   &
848               &                          ' pairing    : '    , TRIM( sdf(jf)%vcomp      ),   &
[1132]849               &                          ' data type: '      ,       sdf(jf)%cltype
[2528]850            call flush(numout)
[1132]851         END DO
852      ENDIF
853     
854   END SUBROUTINE fld_fill
855
856
[1275]857   SUBROUTINE wgt_list( sd, kwgt )
858      !!---------------------------------------------------------------------
859      !!                    ***  ROUTINE wgt_list  ***
860      !!
861      !! ** Purpose :   search array of WGTs and find a weights file
862      !!                entry, or return a new one adding it to the end
863      !!                if it is a new entry, the weights data is read in and
864      !!                restructured (fld_weight)
865      !!----------------------------------------------------------------------
[2715]866      TYPE( FLD ), INTENT(in   ) ::   sd        ! field with name of weights file
867      INTEGER    , INTENT(inout) ::   kwgt      ! index of weights
[1275]868      !!
[2715]869      INTEGER ::   kw, nestid   ! local integer
870      LOGICAL ::   found        ! local logical
[1275]871      !!----------------------------------------------------------------------
872      !
873      !! search down linked list
874      !! weights filename is either present or we hit the end of the list
875      found = .FALSE.
876
877      !! because agrif nest part of filenames are now added in iom_open
878      !! to distinguish between weights files on the different grids, need to track
879      !! nest number explicitly
880      nestid = 0
881#if defined key_agrif
882      nestid = Agrif_Fixed()
883#endif
884      DO kw = 1, nxt_wgt-1
885         IF( TRIM(ref_wgts(kw)%wgtname) == TRIM(sd%wgtname) .AND. &
886             ref_wgts(kw)%nestid == nestid) THEN
887            kwgt = kw
888            found = .TRUE.
889            EXIT
890         ENDIF
891      END DO
892      IF( .NOT.found ) THEN
893         kwgt = nxt_wgt
894         CALL fld_weight( sd )
895      ENDIF
[2715]896      !
[1275]897   END SUBROUTINE wgt_list
898
[2715]899
[1275]900   SUBROUTINE wgt_print( )
901      !!---------------------------------------------------------------------
902      !!                    ***  ROUTINE wgt_print  ***
903      !!
904      !! ** Purpose :   print the list of known weights
905      !!----------------------------------------------------------------------
[2715]906      INTEGER ::   kw   !
[1275]907      !!----------------------------------------------------------------------
908      !
909      DO kw = 1, nxt_wgt-1
910         WRITE(numout,*) 'weight file:  ',TRIM(ref_wgts(kw)%wgtname)
911         WRITE(numout,*) '      ddims:  ',ref_wgts(kw)%ddims(1),ref_wgts(kw)%ddims(2)
912         WRITE(numout,*) '     numwgt:  ',ref_wgts(kw)%numwgt
913         WRITE(numout,*) '     jpiwgt:  ',ref_wgts(kw)%jpiwgt
914         WRITE(numout,*) '     jpjwgt:  ',ref_wgts(kw)%jpjwgt
915         WRITE(numout,*) '    botleft:  ',ref_wgts(kw)%botleft
916         WRITE(numout,*) '   topright:  ',ref_wgts(kw)%topright
917         IF( ref_wgts(kw)%cyclic ) THEN
918            WRITE(numout,*) '       cyclical'
[2528]919            IF( ref_wgts(kw)%overlap > 0 ) WRITE(numout,*) '              with overlap of ', ref_wgts(kw)%overlap
[1275]920         ELSE
921            WRITE(numout,*) '       not cyclical'
922         ENDIF
923         IF( ASSOCIATED(ref_wgts(kw)%data_wgt) )  WRITE(numout,*) '       allocated'
924      END DO
[2715]925      !
[1275]926   END SUBROUTINE wgt_print
927
[2715]928
[1275]929   SUBROUTINE fld_weight( sd )
930      !!---------------------------------------------------------------------
931      !!                    ***  ROUTINE fld_weight  ***
932      !!
933      !! ** Purpose :   create a new WGT structure and fill in data from 
934      !!                file, restructuring as required
935      !!----------------------------------------------------------------------
[2715]936      TYPE( FLD ), INTENT(in) ::   sd   ! field with name of weights file
937      !!
[3294]938      INTEGER                           ::   jn            ! dummy loop indices
939      INTEGER                           ::   inum          ! temporary logical unit
940      INTEGER                           ::   id            ! temporary variable id
941      INTEGER                           ::   ipk           ! temporary vertical dimension
942      CHARACTER (len=5)                 ::   aname
943      INTEGER , DIMENSION(3)            ::   ddims
944      INTEGER , POINTER, DIMENSION(:,:) ::   data_src
945      REAL(wp), POINTER, DIMENSION(:,:) ::   data_tmp
946      LOGICAL                           ::   cyclical
947      INTEGER                           ::   zwrap      ! local integer
[1275]948      !!----------------------------------------------------------------------
949      !
[3294]950      CALL wrk_alloc( jpi,jpj, data_src )   ! integer
951      CALL wrk_alloc( jpi,jpj, data_tmp )
[2715]952      !
[1275]953      IF( nxt_wgt > tot_wgts ) THEN
[2777]954        CALL ctl_stop("fld_weight: weights array size exceeded, increase tot_wgts")
[1275]955      ENDIF
956      !
957      !! new weights file entry, add in extra information
958      !! a weights file represents a 2D grid of a certain shape, so we assume that the current
959      !! input data file is representative of all other files to be opened and processed with the
960      !! current weights file
961
962      !! open input data file (non-model grid)
[1319]963      CALL iom_open( sd%clname, inum, ldiof =  LEN(TRIM(sd%wgtname)) > 0 )
[1275]964
965      !! get dimensions
966      id = iom_varid( inum, sd%clvar, ddims )
967
[2528]968      !! close it
969      CALL iom_close( inum )
[1275]970
[2528]971      !! now open the weights file
[1275]972
[2528]973      CALL iom_open ( sd%wgtname, inum )   ! interpolation weights
974      IF ( inum > 0 ) THEN
[1275]975
[2528]976         !! determine whether we have an east-west cyclic grid
977         !! from global attribute called "ew_wrap" in the weights file
978         !! note that if not found, iom_getatt returns -999 and cyclic with no overlap is assumed
979         !! since this is the most common forcing configuration
[1275]980
[2528]981         CALL iom_getatt(inum, 'ew_wrap', zwrap)
982         IF( zwrap >= 0 ) THEN
[1275]983            cyclical = .TRUE.
[2528]984         ELSE IF( zwrap == -999 ) THEN
[1275]985            cyclical = .TRUE.
[2528]986            zwrap = 0
987         ELSE
988            cyclical = .FALSE.
[1275]989         ENDIF
990
991         ref_wgts(nxt_wgt)%ddims(1) = ddims(1)
992         ref_wgts(nxt_wgt)%ddims(2) = ddims(2)
993         ref_wgts(nxt_wgt)%wgtname = sd%wgtname
[2528]994         ref_wgts(nxt_wgt)%overlap = zwrap
995         ref_wgts(nxt_wgt)%cyclic = cyclical
[1275]996         ref_wgts(nxt_wgt)%nestid = 0
997#if defined key_agrif
998         ref_wgts(nxt_wgt)%nestid = Agrif_Fixed()
999#endif
1000         !! weights file is stored as a set of weights (wgt01->wgt04 or wgt01->wgt16)
1001         !! for each weight wgtNN there is an integer array srcNN which gives the point in
1002         !! the input data grid which is to be multiplied by the weight
1003         !! they are both arrays on the model grid so the result of the multiplication is
1004         !! added into an output array on the model grid as a running sum
1005
1006         !! two possible cases: bilinear (4 weights) or bicubic (16 weights)
1007         id = iom_varid(inum, 'src05', ldstop=.FALSE.)
1008         IF( id <= 0) THEN
1009            ref_wgts(nxt_wgt)%numwgt = 4
1010         ELSE
1011            ref_wgts(nxt_wgt)%numwgt = 16
1012         ENDIF
1013
1014         ALLOCATE( ref_wgts(nxt_wgt)%data_jpi(jpi,jpj,4) )
1015         ALLOCATE( ref_wgts(nxt_wgt)%data_jpj(jpi,jpj,4) )
1016         ALLOCATE( ref_wgts(nxt_wgt)%data_wgt(jpi,jpj,ref_wgts(nxt_wgt)%numwgt) )
1017
1018         DO jn = 1,4
1019            aname = ' '
1020            WRITE(aname,'(a3,i2.2)') 'src',jn
1021            data_tmp(:,:) = 0
[1955]1022            CALL iom_get ( inum, jpdom_data, aname, data_tmp(:,:) )
[1275]1023            data_src(:,:) = INT(data_tmp(:,:))
1024            ref_wgts(nxt_wgt)%data_jpj(:,:,jn) = 1 + (data_src(:,:)-1) / ref_wgts(nxt_wgt)%ddims(1)
1025            ref_wgts(nxt_wgt)%data_jpi(:,:,jn) = data_src(:,:) - ref_wgts(nxt_wgt)%ddims(1)*(ref_wgts(nxt_wgt)%data_jpj(:,:,jn)-1)
1026         END DO
1027
1028         DO jn = 1, ref_wgts(nxt_wgt)%numwgt
1029            aname = ' '
1030            WRITE(aname,'(a3,i2.2)') 'wgt',jn
[1955]1031            ref_wgts(nxt_wgt)%data_wgt(:,:,jn) = 0.0
1032            CALL iom_get ( inum, jpdom_data, aname, ref_wgts(nxt_wgt)%data_wgt(:,:,jn) )
[1275]1033         END DO
1034         CALL iom_close (inum)
1035 
1036         ! find min and max indices in grid
[1955]1037         ref_wgts(nxt_wgt)%botleft(1) = MINVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1038         ref_wgts(nxt_wgt)%botleft(2) = MINVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1039         ref_wgts(nxt_wgt)%topright(1) = MAXVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1040         ref_wgts(nxt_wgt)%topright(2) = MAXVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
[1275]1041
1042         ! and therefore dimensions of the input box
1043         ref_wgts(nxt_wgt)%jpiwgt = ref_wgts(nxt_wgt)%topright(1) - ref_wgts(nxt_wgt)%botleft(1) + 1
1044         ref_wgts(nxt_wgt)%jpjwgt = ref_wgts(nxt_wgt)%topright(2) - ref_wgts(nxt_wgt)%botleft(2) + 1
1045
1046         ! shift indexing of source grid
1047         ref_wgts(nxt_wgt)%data_jpi(:,:,:) = ref_wgts(nxt_wgt)%data_jpi(:,:,:) - ref_wgts(nxt_wgt)%botleft(1) + 1
1048         ref_wgts(nxt_wgt)%data_jpj(:,:,:) = ref_wgts(nxt_wgt)%data_jpj(:,:,:) - ref_wgts(nxt_wgt)%botleft(2) + 1
1049
1050         ! create input grid, give it a halo to allow gradient calculations
[1702]1051         ! SA: +3 stencil is a patch to avoid out-of-bound computation in some configuration.
1052         ! a more robust solution will be given in next release
[2528]1053         ipk =  SIZE(sd%fnow, 3)
1054         ALLOCATE( ref_wgts(nxt_wgt)%fly_dta(ref_wgts(nxt_wgt)%jpiwgt+3, ref_wgts(nxt_wgt)%jpjwgt+3 ,ipk) )
1055         IF( ref_wgts(nxt_wgt)%cyclic ) ALLOCATE( ref_wgts(nxt_wgt)%col(1,ref_wgts(nxt_wgt)%jpjwgt+3,ipk) )
[1275]1056
1057         nxt_wgt = nxt_wgt + 1
1058
1059      ELSE
1060         CALL ctl_stop( '    fld_weight : unable to read the file ' )
1061      ENDIF
1062
[3294]1063      CALL wrk_dealloc( jpi,jpj, data_src )   ! integer
1064      CALL wrk_dealloc( jpi,jpj, data_tmp )
[2715]1065      !
[1275]1066   END SUBROUTINE fld_weight
1067
[2715]1068
1069   SUBROUTINE fld_interp( num, clvar, kw, kk, dta, nrec )
[1275]1070      !!---------------------------------------------------------------------
1071      !!                    ***  ROUTINE fld_interp  ***
1072      !!
1073      !! ** Purpose :   apply weights to input gridded data to create data
1074      !!                on model grid
1075      !!----------------------------------------------------------------------
[2715]1076      INTEGER                   , INTENT(in   ) ::   num     ! stream number
1077      CHARACTER(LEN=*)          , INTENT(in   ) ::   clvar   ! variable name
1078      INTEGER                   , INTENT(in   ) ::   kw      ! weights number
1079      INTEGER                   , INTENT(in   ) ::   kk      ! vertical dimension of kk
1080      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   dta     ! output field on model grid
1081      INTEGER                   , INTENT(in   ) ::   nrec    ! record number to read (ie time slice)
[1275]1082      !!
[2715]1083      INTEGER, DIMENSION(3) ::   rec1,recn   ! temporary arrays for start and length
1084      INTEGER ::  jk, jn, jm           ! loop counters
1085      INTEGER ::  ni, nj               ! lengths
1086      INTEGER ::  jpimin,jpiwid        ! temporary indices
1087      INTEGER ::  jpjmin,jpjwid        ! temporary indices
1088      INTEGER ::  jpi1,jpi2,jpj1,jpj2  ! temporary indices
[1275]1089      !!----------------------------------------------------------------------
1090      !
1091      !! for weighted interpolation we have weights at four corners of a box surrounding
1092      !! a model grid point, each weight is multiplied by a grid value (bilinear case)
1093      !! or by a grid value and gradients at the corner point (bicubic case)
1094      !! so we need to have a 4 by 4 subgrid surrounding each model point to cover both cases
1095
[2528]1096      !! sub grid from non-model input grid which encloses all grid points in this nemo process
[1275]1097      jpimin = ref_wgts(kw)%botleft(1)
1098      jpjmin = ref_wgts(kw)%botleft(2)
1099      jpiwid = ref_wgts(kw)%jpiwgt
1100      jpjwid = ref_wgts(kw)%jpjwgt
1101
[2528]1102      !! when reading in, expand this sub-grid by one halo point all the way round for calculating gradients
[1275]1103      rec1(1) = MAX( jpimin-1, 1 )
1104      rec1(2) = MAX( jpjmin-1, 1 )
[2528]1105      rec1(3) = 1
[1275]1106      recn(1) = MIN( jpiwid+2, ref_wgts(kw)%ddims(1)-rec1(1)+1 )
1107      recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
[2528]1108      recn(3) = kk
[1275]1109
[2528]1110      !! where we need to put it in the non-nemo grid fly_dta
1111      !! note that jpi1 and jpj1 only differ from 1 when jpimin and jpjmin are 1
1112      !! (ie at the extreme west or south of the whole input grid) and similarly for jpi2 and jpj2
[1275]1113      jpi1 = 2 + rec1(1) - jpimin
1114      jpj1 = 2 + rec1(2) - jpjmin
1115      jpi2 = jpi1 + recn(1) - 1
1116      jpj2 = jpj1 + recn(2) - 1
1117
[2528]1118      ref_wgts(kw)%fly_dta(:,:,:) = 0.0
1119      SELECT CASE( SIZE(ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:),3) )
1120      CASE(1)
1121           CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,1), nrec, rec1, recn)
1122      CASE DEFAULT
1123           CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:), nrec, rec1, recn)
1124      END SELECT 
[1275]1125
1126      !! first four weights common to both bilinear and bicubic
[2528]1127      !! data_jpi, data_jpj have already been shifted to (1,1) corresponding to botleft
[1275]1128      !! note that we have to offset by 1 into fly_dta array because of halo
[2528]1129      dta(:,:,:) = 0.0
[1275]1130      DO jk = 1,4
1131        DO jn = 1, jpj
1132          DO jm = 1,jpi
1133            ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1134            nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
[2528]1135            dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk) * ref_wgts(kw)%fly_dta(ni+1,nj+1,:)
[1275]1136          END DO
1137        END DO
1138      END DO
1139
1140      IF (ref_wgts(kw)%numwgt .EQ. 16) THEN
1141
1142        !! fix up halo points that we couldnt read from file
1143        IF( jpi1 == 2 ) THEN
[2528]1144           ref_wgts(kw)%fly_dta(jpi1-1,:,:) = ref_wgts(kw)%fly_dta(jpi1,:,:)
[1275]1145        ENDIF
1146        IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
[2528]1147           ref_wgts(kw)%fly_dta(jpi2+1,:,:) = ref_wgts(kw)%fly_dta(jpi2,:,:)
[1275]1148        ENDIF
1149        IF( jpj1 == 2 ) THEN
[2528]1150           ref_wgts(kw)%fly_dta(:,jpj1-1,:) = ref_wgts(kw)%fly_dta(:,jpj1,:)
[1275]1151        ENDIF
1152        IF( jpj2 + jpjmin - 1 == ref_wgts(kw)%ddims(2)+1 .AND. jpj2 .lt. jpjwid+2 ) THEN
[2528]1153           ref_wgts(kw)%fly_dta(:,jpj2+1,:) = 2.0*ref_wgts(kw)%fly_dta(:,jpj2,:) - ref_wgts(kw)%fly_dta(:,jpj2-1,:)
[1275]1154        ENDIF
1155
1156        !! if data grid is cyclic we can do better on east-west edges
1157        !! but have to allow for whether first and last columns are coincident
1158        IF( ref_wgts(kw)%cyclic ) THEN
1159           rec1(2) = MAX( jpjmin-1, 1 )
[2528]1160           recn(1) = 1
[1275]1161           recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
1162           jpj1 = 2 + rec1(2) - jpjmin
1163           jpj2 = jpj1 + recn(2) - 1
1164           IF( jpi1 == 2 ) THEN
[2528]1165              rec1(1) = ref_wgts(kw)%ddims(1) - ref_wgts(kw)%overlap
1166              SELECT CASE( SIZE( ref_wgts(kw)%col(:,jpj1:jpj2,:),3) )
1167              CASE(1)
1168                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,1), nrec, rec1, recn)
1169              CASE DEFAULT
1170                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, rec1, recn)
1171              END SELECT     
1172              ref_wgts(kw)%fly_dta(jpi1-1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
[1275]1173           ENDIF
1174           IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
[2528]1175              rec1(1) = 1 + ref_wgts(kw)%overlap
1176              SELECT CASE( SIZE( ref_wgts(kw)%col(:,jpj1:jpj2,:),3) )
1177              CASE(1)
1178                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,1), nrec, rec1, recn)
1179              CASE DEFAULT
1180                   CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, rec1, recn)
1181              END SELECT
1182              ref_wgts(kw)%fly_dta(jpi2+1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
[1275]1183           ENDIF
1184        ENDIF
1185
1186        ! gradient in the i direction
1187        DO jk = 1,4
1188          DO jn = 1, jpj
1189            DO jm = 1,jpi
1190              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1191              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
[2528]1192              dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+4) * 0.5 *         &
1193                               (ref_wgts(kw)%fly_dta(ni+2,nj+1,:) - ref_wgts(kw)%fly_dta(ni,nj+1,:))
[1275]1194            END DO
1195          END DO
1196        END DO
1197
1198        ! gradient in the j direction
1199        DO jk = 1,4
1200          DO jn = 1, jpj
1201            DO jm = 1,jpi
1202              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1203              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
[2528]1204              dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+8) * 0.5 *         &
1205                               (ref_wgts(kw)%fly_dta(ni+1,nj+2,:) - ref_wgts(kw)%fly_dta(ni+1,nj,:))
[1275]1206            END DO
1207          END DO
1208        END DO
1209
[2715]1210         ! gradient in the ij direction
1211         DO jk = 1,4
1212            DO jn = 1, jpj
1213               DO jm = 1,jpi
1214                  ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
1215                  nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
1216                  dta(jm,jn,:) = dta(jm,jn,:) + ref_wgts(kw)%data_wgt(jm,jn,jk+12) * 0.25 * ( &
[2528]1217                               (ref_wgts(kw)%fly_dta(ni+2,nj+2,:) - ref_wgts(kw)%fly_dta(ni  ,nj+2,:)) -   &
1218                               (ref_wgts(kw)%fly_dta(ni+2,nj  ,:) - ref_wgts(kw)%fly_dta(ni  ,nj  ,:)))
[2715]1219               END DO
[1275]1220            END DO
[2715]1221         END DO
1222         !
[1275]1223      END IF
[2715]1224      !
[1275]1225   END SUBROUTINE fld_interp
[2528]1226
1227
1228   FUNCTION ksec_week( cdday )
1229      !!---------------------------------------------------------------------
1230      !!                    ***  FUNCTION kshift_week ***
1231      !!
1232      !! ** Purpose : 
1233      !!---------------------------------------------------------------------
1234      CHARACTER(len=*), INTENT(in)   ::   cdday   !3 first letters of the first day of the weekly file
1235      !!
1236      INTEGER                        ::   ksec_week  ! output variable
1237      INTEGER                        ::   ijul       !temp variable
1238      INTEGER                        ::   ishift     !temp variable
1239      CHARACTER(len=3),DIMENSION(7)  ::   cl_week 
1240      !!----------------------------------------------------------------------
1241      cl_week = (/"sun","sat","fri","thu","wed","tue","mon"/)
1242      DO ijul = 1, 7
1243         IF( cl_week(ijul) == TRIM(cdday) ) EXIT
[2715]1244      END DO
[2528]1245      IF( ijul .GT. 7 )   CALL ctl_stop( 'ksec_week: wrong day for sdjf%cltype(6:8): '//TRIM(cdday) )
1246      !
1247      ishift = ijul * NINT(rday)
1248      !
1249      ksec_week = nsec_week + ishift
1250      ksec_week = MOD( ksec_week, 7*NINT(rday) )
1251      !
1252   END FUNCTION ksec_week
1253
[2715]1254   !!======================================================================
[888]1255END MODULE fldread
Note: See TracBrowser for help on using the repository browser.