New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
fldread.F90 in NEMO/branches/UKMO/dev_r12866_HPC-02_Daley_Tiling_trial_extra_halo/src/OCE/SBC – NEMO

source: NEMO/branches/UKMO/dev_r12866_HPC-02_Daley_Tiling_trial_extra_halo/src/OCE/SBC/fldread.F90 @ 12962

Last change on this file since 12962 was 12962, checked in by hadcv, 4 years ago

Update with [12960] fixes

  • Property svn:keywords set to Id
File size: 88.2 KB
Line 
1MODULE fldread
2   !!======================================================================
3   !!                       ***  MODULE  fldread  ***
4   !! Ocean forcing:  read input field for surface boundary condition
5   !!=====================================================================
6   !! History :  2.0  !  2006-06  (S. Masson, G. Madec)  Original code
7   !!            3.0  !  2008-05  (S. Alderson)  Modified for Interpolation in memory from input grid to model grid
8   !!            3.4  !  2013-10  (D. Delrosso, P. Oddo)  suppression of land point prior to interpolation
9   !!                 !  12-2015  (J. Harle) Adding BDY on-the-fly interpolation
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   fld_read      : read input fields used for the computation of the surface boundary condition
14   !!   fld_init      : initialization of field read
15   !!   fld_def       : define the record(s) of the file and its name
16   !!   fld_get       : read the data
17   !!   fld_map       : read global data from file and map onto local data using a general mapping (use for open boundaries)
18   !!   fld_rot       : rotate the vector fields onto the local grid direction
19   !!   fld_clopn     : close/open the files
20   !!   fld_fill      : fill the data structure with the associated information read in namelist
21   !!   wgt_list      : manage the weights used for interpolation
22   !!   wgt_print     : print the list of known weights
23   !!   fld_weight    : create a WGT structure and fill in data from file, restructuring as required
24   !!   apply_seaoverland : fill land with ocean values
25   !!   seaoverland   : create shifted matrices for seaoverland application
26   !!   fld_interp    : apply weights to input gridded data to create data on model grid
27   !!   fld_filename  : define the filename according to a given date
28   !!   ksec_week     : function returning seconds between 00h of the beginning of the week and half of the current time step
29   !!----------------------------------------------------------------------
30   USE oce            ! ocean dynamics and tracers
31   USE dom_oce        ! ocean space and time domain
32   USE phycst         ! physical constant
33   USE sbc_oce        ! surface boundary conditions : fields
34   USE geo2ocean      ! for vector rotation on to model grid
35   !
36   USE in_out_manager ! I/O manager
37   USE iom            ! I/O manager library
38   USE ioipsl  , ONLY : ymds2ju, ju2ymds   ! for calendar
39   USE lib_mpp        ! MPP library
40   USE lbclnk         ! ocean lateral boundary conditions (C1D case)
41   
42   IMPLICIT NONE
43   PRIVATE   
44 
45   PUBLIC   fld_map    ! routine called by tides_init
46   PUBLIC   fld_read, fld_fill   ! called by sbc... modules
47   PUBLIC   fld_def
48
49   TYPE, PUBLIC ::   FLD_N      !: Namelist field informations
50      CHARACTER(len = 256) ::   clname      ! generic name of the NetCDF flux file
51      REAL(wp)             ::   freqh       ! frequency of each flux file
52      CHARACTER(len = 34)  ::   clvar       ! generic name of the variable in the NetCDF flux file
53      LOGICAL              ::   ln_tint     ! time interpolation or not (T/F)
54      LOGICAL              ::   ln_clim     ! climatology or not (T/F)
55      CHARACTER(len = 8)   ::   clftyp      ! type of data file 'daily', 'monthly' or yearly'
56      CHARACTER(len = 256) ::   wname       ! generic name of a NetCDF weights file to be used, blank if not
57      CHARACTER(len = 34)  ::   vcomp       ! symbolic component name if a vector that needs rotation
58      !                                     ! a string starting with "U" or "V" for each component   
59      !                                     ! chars 2 onwards identify which components go together 
60      CHARACTER(len = 34)  ::   lname       ! generic name of a NetCDF land/sea mask file to be used, blank if not
61      !                                     ! 0=sea 1=land
62   END TYPE FLD_N
63
64   TYPE, PUBLIC ::   FLD        !: Input field related variables
65      CHARACTER(len = 256)            ::   clrootname   ! generic name of the NetCDF file
66      CHARACTER(len = 256)            ::   clname       ! current name of the NetCDF file
67      REAL(wp)                        ::   freqh        ! frequency of each flux file
68      CHARACTER(len = 34)             ::   clvar        ! generic name of the variable in the NetCDF flux file
69      LOGICAL                         ::   ln_tint      ! time interpolation or not (T/F)
70      LOGICAL                         ::   ln_clim      ! climatology or not (T/F)
71      CHARACTER(len = 8)              ::   clftyp       ! type of data file 'daily', 'monthly' or yearly'
72      CHARACTER(len = 1)              ::   cltype       ! nature of grid-points: T, U, V...
73      REAL(wp)                        ::   zsgn         ! -1. the sign change across the north fold, =  1. otherwise
74      INTEGER                         ::   num          ! iom id of the jpfld files to be read
75      INTEGER , DIMENSION(2,2)        ::   nrec         ! before/after record (1: index, 2: second since Jan. 1st 00h of yr nit000)
76      INTEGER                         ::   nbb          ! index of before values
77      INTEGER                         ::   naa          ! index of after  values
78      INTEGER , ALLOCATABLE, DIMENSION(:) ::   nrecsec   !
79      REAL(wp), POINTER, DIMENSION(:,:,:  ) ::   fnow   ! input fields interpolated to now time step
80      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   fdta   ! 2 consecutive record of input fields
81      CHARACTER(len = 256)            ::   wgtname      ! current name of the NetCDF weight file acting as a key
82      !                                                 ! into the WGTLIST structure
83      CHARACTER(len = 34)             ::   vcomp        ! symbolic name for a vector component that needs rotation
84      LOGICAL, DIMENSION(2)           ::   rotn         ! flag to indicate whether before/after field has been rotated
85      INTEGER                         ::   nreclast     ! last record to be read in the current file
86      CHARACTER(len = 256)            ::   lsmname      ! current name of the NetCDF mask file acting as a key
87      !                                                 !
88      !                                                 ! Variables related to BDY
89      INTEGER                         ::   igrd         !   grid type for bdy data
90      INTEGER                         ::   ibdy         !   bdy set id number
91      INTEGER, POINTER, DIMENSION(:)  ::   imap         !   Array of integer pointers to 1D arrays
92      LOGICAL                         ::   ltotvel      !   total velocity or not (T/F)
93      LOGICAL                         ::   lzint        !   T if it requires a vertical interpolation
94   END TYPE FLD
95
96!$AGRIF_DO_NOT_TREAT
97
98   !! keep list of all weights variables so they're only read in once
99   !! need to add AGRIF directives not to process this structure
100   !! also need to force wgtname to include AGRIF nest number
101   TYPE         ::   WGT        !: Input weights related variables
102      CHARACTER(len = 256)                    ::   wgtname      ! current name of the NetCDF weight file
103      INTEGER , DIMENSION(2)                  ::   ddims        ! shape of input grid
104      INTEGER , DIMENSION(2)                  ::   botleft      ! top left corner of box in input grid containing
105      !                                                         ! current processor grid
106      INTEGER , DIMENSION(2)                  ::   topright     ! top right corner of box
107      INTEGER                                 ::   jpiwgt       ! width of box on input grid
108      INTEGER                                 ::   jpjwgt       ! height of box on input grid
109      INTEGER                                 ::   numwgt       ! number of weights (4=bilinear, 16=bicubic)
110      INTEGER                                 ::   nestid       ! for agrif, keep track of nest we're in
111      INTEGER                                 ::   overlap      ! =0 when cyclic grid has no overlapping EW columns
112      !                                                         ! =>1 when they have one or more overlapping columns     
113      !                                                         ! =-1 not cyclic
114      LOGICAL                                 ::   cyclic       ! east-west cyclic or not
115      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpi     ! array of source integers
116      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpj     ! array of source integers
117      REAL(wp), DIMENSION(:,:,:), POINTER     ::   data_wgt     ! array of weights on model grid
118      REAL(wp), DIMENSION(:,:,:), POINTER     ::   fly_dta      ! array of values on input grid
119      REAL(wp), DIMENSION(:,:,:), POINTER     ::   col          ! temporary array for reading in columns
120   END TYPE WGT
121
122   INTEGER,     PARAMETER             ::   tot_wgts = 20
123   TYPE( WGT ), DIMENSION(tot_wgts)   ::   ref_wgts     ! array of wgts
124   INTEGER                            ::   nxt_wgt = 1  ! point to next available space in ref_wgts array
125   INTEGER                            ::   nflag = 0
126   REAL(wp), PARAMETER                ::   undeff_lsm = -999.00_wp
127
128!$AGRIF_END_DO_NOT_TREAT
129
130   !! * Substitutions
131#  include "do_loop_substitute.h90"
132   !!----------------------------------------------------------------------
133   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
134   !! $Id$
135   !! Software governed by the CeCILL license (see ./LICENSE)
136   !!----------------------------------------------------------------------
137CONTAINS
138
139   SUBROUTINE fld_read( kt, kn_fsbc, sd, kit, pt_offset, Kmm )
140      !!---------------------------------------------------------------------
141      !!                    ***  ROUTINE fld_read  ***
142      !!                   
143      !! ** Purpose :   provide at each time step the surface ocean fluxes
144      !!                (momentum, heat, freshwater and runoff)
145      !!
146      !! ** Method  :   READ each input fields in NetCDF files using IOM
147      !!      and intepolate it to the model time-step.
148      !!         Several assumptions are made on the input file:
149      !!      blahblahblah....
150      !!----------------------------------------------------------------------
151      INTEGER  , INTENT(in   )               ::   kt        ! ocean time step
152      INTEGER  , INTENT(in   )               ::   kn_fsbc   ! sbc computation period (in time step)
153      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
154      INTEGER  , INTENT(in   ), OPTIONAL     ::   kit       ! subcycle timestep for timesplitting option
155      REAL(wp) , INTENT(in   ), OPTIONAL     ::   pt_offset ! provide fields at time other than "now"
156      INTEGER  , INTENT(in   ), OPTIONAL     ::   Kmm       ! ocean time level index
157      !!
158      INTEGER  ::   imf          ! size of the structure sd
159      INTEGER  ::   jf           ! dummy indices
160      INTEGER  ::   isecsbc      ! number of seconds between Jan. 1st 00h of nit000 year and the middle of sbc time step
161      INTEGER  ::   ibb, iaa     ! shorter name for sd(jf)%nbb and sd(jf)%naa
162      LOGICAL  ::   ll_firstcall ! true if this is the first call to fld_read for this set of fields
163      REAL(wp) ::   zt_offset    ! local time offset variable
164      REAL(wp) ::   ztinta       ! ratio applied to after  records when doing time interpolation
165      REAL(wp) ::   ztintb       ! ratio applied to before records when doing time interpolation
166      CHARACTER(LEN=1000) ::   clfmt  ! write format
167      !!---------------------------------------------------------------------
168      ll_firstcall = kt == nit000
169      IF( PRESENT(kit) )   ll_firstcall = ll_firstcall .and. kit == 1
170
171      IF( nn_components == jp_iam_sas ) THEN   ;   zt_offset = REAL( nn_fsbc, wp )
172      ELSE                                     ;   zt_offset = 0.
173      ENDIF
174      IF( PRESENT(pt_offset) )   zt_offset = pt_offset
175
176      ! Note that all varibles starting by nsec_* are shifted time by +1/2 time step to be centrered
177      IF( PRESENT(kit) ) THEN   ! ignore kn_fsbc in this case
178         isecsbc = nsec_year + nsec1jan000 + NINT( (     REAL(      kit,wp) + zt_offset ) * rn_Dt / REAL(nn_e,wp) )
179      ELSE                      ! middle of sbc time step
180         ! note: we use kn_fsbc-1 because nsec_year is defined at the middle of the current time step
181         isecsbc = nsec_year + nsec1jan000 + NINT( ( 0.5*REAL(kn_fsbc-1,wp) + zt_offset ) * rn_Dt )
182      ENDIF
183      imf = SIZE( sd )
184      !
185      IF( ll_firstcall ) THEN                      ! initialization
186         DO jf = 1, imf 
187            IF( TRIM(sd(jf)%clrootname) == 'NOT USED' )   CYCLE
188            CALL fld_init( isecsbc, sd(jf) )       ! read each before field (put them in after as they will be swapped)
189         END DO
190         IF( lwp ) CALL wgt_print()                ! control print
191      ENDIF
192      !                                            ! ====================================== !
193      IF( MOD( kt-1, kn_fsbc ) == 0 ) THEN         ! update field at each kn_fsbc time-step !
194         !                                         ! ====================================== !
195         !
196         DO jf = 1, imf                            ! ---   loop over field   --- !
197            !
198            IF( TRIM(sd(jf)%clrootname) == 'NOT USED' )   CYCLE
199            CALL fld_update( isecsbc, sd(jf), Kmm )
200            !
201         END DO                                    ! --- end loop over field --- !
202
203         CALL fld_rot( kt, sd )                    ! rotate vector before/now/after fields if needed
204
205         DO jf = 1, imf                            ! ---   loop over field   --- !
206            !
207            IF( TRIM(sd(jf)%clrootname) == 'NOT USED' )   CYCLE
208            !
209            ibb = sd(jf)%nbb   ;   iaa = sd(jf)%naa
210            !
211            IF( sd(jf)%ln_tint ) THEN              ! temporal interpolation
212               IF(lwp .AND. kt - nit000 <= 100 ) THEN
213                  clfmt = "('   fld_read: var ', a, ' kt = ', i8, ' (', f9.4,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
214                     &    "', records b/a: ', i6.4, '/', i6.4, ' (days ', f9.4,'/', f9.4, ')')"
215                  WRITE(numout, clfmt)  TRIM( sd(jf)%clvar ), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,   &           
216                     & sd(jf)%nrec(1,ibb), sd(jf)%nrec(1,iaa), REAL(sd(jf)%nrec(2,ibb),wp)/rday, REAL(sd(jf)%nrec(2,iaa),wp)/rday
217                  WRITE(numout, *) '      zt_offset is : ',zt_offset
218               ENDIF
219               ! temporal interpolation weights
220               ztinta =  REAL( isecsbc - sd(jf)%nrec(2,ibb), wp ) / REAL( sd(jf)%nrec(2,iaa) - sd(jf)%nrec(2,ibb), wp )
221               ztintb =  1. - ztinta
222               sd(jf)%fnow(:,:,:) = ztintb * sd(jf)%fdta(:,:,:,ibb) + ztinta * sd(jf)%fdta(:,:,:,iaa)
223            ELSE   ! nothing to do...
224               IF(lwp .AND. kt - nit000 <= 100 ) THEN
225                  clfmt = "('   fld_read: var ', a, ' kt = ', i8,' (', f9.4,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
226                     &    "', record: ', i6.4, ' (days ', f9.4, ' <-> ', f9.4, ')')"
227                  WRITE(numout, clfmt) TRIM(sd(jf)%clvar), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,    &
228                     &                 sd(jf)%nrec(1,iaa), REAL(sd(jf)%nrec(2,ibb),wp)/rday, REAL(sd(jf)%nrec(2,iaa),wp)/rday
229               ENDIF
230            ENDIF
231            !
232            IF( kt == nitend - kn_fsbc + 1 )   CALL iom_close( sd(jf)%num )   ! Close the input files
233
234         END DO                                    ! --- end loop over field --- !
235         !
236      ENDIF
237      !
238   END SUBROUTINE fld_read
239
240
241   SUBROUTINE fld_init( ksecsbc, sdjf )
242      !!---------------------------------------------------------------------
243      !!                    ***  ROUTINE fld_init  ***
244      !!
245      !! ** Purpose :  - first call(s) to fld_def to define before values
246      !!               - open file
247      !!----------------------------------------------------------------------
248      INTEGER  , INTENT(in   ) ::   ksecsbc   !
249      TYPE(FLD), INTENT(inout) ::   sdjf         ! input field related variables
250      !!---------------------------------------------------------------------
251      !
252      IF( nflag == 0 )   nflag = -( HUGE(0) - 10 )
253      !
254      CALL fld_def( sdjf )
255      IF( sdjf%ln_tint .AND. ksecsbc < sdjf%nrecsec(1) )   CALL fld_def( sdjf, ldprev = .TRUE. )
256      !
257      CALL fld_clopn( sdjf )
258      sdjf%nrec(:,sdjf%naa) = (/ 1, nflag /)  ! default definition to force flp_update to read the file.
259      !
260   END SUBROUTINE fld_init
261
262
263   SUBROUTINE fld_update( ksecsbc, sdjf, Kmm )
264      !!---------------------------------------------------------------------
265      !!                    ***  ROUTINE fld_update  ***
266      !!
267      !! ** Purpose : Compute
268      !!              if sdjf%ln_tint = .TRUE.
269      !!                  nrec(:,iaa): record number and its time (nrec(:,ibb) is obtained from nrec(:,iaa) when swapping)
270      !!              if sdjf%ln_tint = .FALSE.
271      !!                  nrec(1,iaa): record number
272      !!                  nrec(2,ibb) and nrec(2,iaa): time of the beginning and end of the record
273      !!----------------------------------------------------------------------
274      INTEGER  ,           INTENT(in   ) ::   ksecsbc   !
275      TYPE(FLD),           INTENT(inout) ::   sdjf      ! input field related variables
276      INTEGER  , OPTIONAL, INTENT(in   ) ::   Kmm    ! ocean time level index
277      !
278      INTEGER  ::   ja           ! end of this record (in seconds)
279      INTEGER  ::   ibb, iaa     ! shorter name for sdjf%nbb and sdjf%naa
280      !!----------------------------------------------------------------------
281      ibb = sdjf%nbb   ;   iaa = sdjf%naa
282      !
283      IF( ksecsbc > sdjf%nrec(2,iaa) ) THEN     ! --> we need to update after data
284       
285         ! find where is the new after record... (it is not necessary sdjf%nrec(1,iaa)+1 )
286         ja = sdjf%nrec(1,iaa)
287         DO WHILE ( ksecsbc >= sdjf%nrecsec(ja) .AND. ja < sdjf%nreclast )   ! Warning: make sure ja <= sdjf%nreclast in this test
288            ja = ja + 1
289         END DO
290         IF( ksecsbc > sdjf%nrecsec(ja) )   ja = ja + 1   ! in case ksecsbc > sdjf%nrecsec(sdjf%nreclast)
291
292         ! if ln_tint and if the new after is not ja+1, we need also to update after data before the swap
293         ! so, after the swap, sdjf%nrec(2,ibb) will still be the closest value located just before ksecsbc
294         IF( sdjf%ln_tint .AND. ( ja > sdjf%nrec(1,iaa) + 1 .OR. sdjf%nrec(2,iaa) == nflag ) ) THEN
295            sdjf%nrec(:,iaa) = (/ ja-1, sdjf%nrecsec(ja-1) /)   ! update nrec(:,iaa) with before information
296            CALL fld_get( sdjf, Kmm )                           ! read after data that will be used as before data
297         ENDIF
298           
299         ! if after is in the next file...
300         IF( ja > sdjf%nreclast ) THEN
301           
302            CALL fld_def( sdjf )
303            IF( ksecsbc > sdjf%nrecsec(sdjf%nreclast) )   CALL fld_def( sdjf, ldnext = .TRUE. )
304            CALL fld_clopn( sdjf )           ! open next file
305           
306            ! find where is after in this new file
307            ja = 1
308            DO WHILE ( ksecsbc > sdjf%nrecsec(ja) .AND. ja < sdjf%nreclast )
309               ja = ja + 1
310            END DO
311            IF( ksecsbc > sdjf%nrecsec(ja) )   ja = ja + 1   ! in case ksecsbc > sdjf%nrecsec(sdjf%nreclast)
312           
313            IF( ja > sdjf%nreclast ) THEN
314               CALL ctl_stop( "STOP", "fld_def: need next-next file? we should not be there... file: "//TRIM(sdjf%clrootname) )
315            ENDIF
316           
317            ! if ln_tint and if after is not the first record, we must (potentially again) update after data before the swap
318            IF( sdjf%ln_tint .AND. ja > 1 ) THEN
319               IF( sdjf%nrecsec(0) /= nflag ) THEN                    ! no trick used: after file is not the current file
320                  sdjf%nrec(:,iaa) = (/ ja-1, sdjf%nrecsec(ja-1) /)   ! update nrec(:,iaa) with before information
321                  CALL fld_get( sdjf, Kmm )                           ! read after data that will be used as before data
322               ENDIF
323            ENDIF
324           
325         ENDIF
326
327         IF( sdjf%ln_tint ) THEN                                ! Swap data
328            sdjf%nbb = sdjf%naa                                 !    swap indices
329            sdjf%naa = 3 - sdjf%naa                             !    = 2(1) if naa == 1(2)
330         ELSE                                                   ! No swap
331            sdjf%nrec(:,ibb) = (/ ja-1, sdjf%nrecsec(ja-1) /)   !    only for print
332         ENDIF
333           
334         ! read new after data
335         sdjf%nrec(:,sdjf%naa) = (/ ja, sdjf%nrecsec(ja) /)     ! update nrec(:,naa) as it is used by fld_get
336         CALL fld_get( sdjf, Kmm )                              ! read after data (with nrec(:,naa) informations)
337       
338      ENDIF
339      !
340   END SUBROUTINE fld_update
341
342
343   SUBROUTINE fld_get( sdjf, Kmm )
344      !!---------------------------------------------------------------------
345      !!                    ***  ROUTINE fld_get  ***
346      !!
347      !! ** Purpose :   read the data
348      !!----------------------------------------------------------------------
349      TYPE(FLD),           INTENT(inout) ::   sdjf   ! input field related variables
350      INTEGER  , OPTIONAL, INTENT(in   ) ::   Kmm    ! ocean time level index
351      !
352      INTEGER ::   ipk      ! number of vertical levels of sdjf%fdta ( 2D: ipk=1 ; 3D: ipk=jpk )
353      INTEGER ::   iaa      ! shorter name for sdjf%naa
354      INTEGER ::   iw       ! index into wgts array
355      INTEGER ::   idvar    ! variable ID
356      INTEGER ::   idmspc   ! number of spatial dimensions
357      LOGICAL ::   lmoor    ! C1D case: point data
358      REAL(wp), DIMENSION(:,:,:), POINTER ::   dta_alias   ! short cut
359      !!---------------------------------------------------------------------
360      iaa = sdjf%naa
361      !
362      IF( sdjf%ln_tint ) THEN   ;   dta_alias => sdjf%fdta(:,:,:,iaa)
363      ELSE                      ;   dta_alias => sdjf%fnow(:,:,:    )
364      ENDIF
365      ipk = SIZE( dta_alias, 3 )
366      !
367      IF( ASSOCIATED(sdjf%imap) ) THEN              ! BDY case
368         CALL fld_map( sdjf%num, sdjf%clvar, dta_alias(:,:,:), sdjf%nrec(1,iaa),   &
369            &          sdjf%imap, sdjf%igrd, sdjf%ibdy, sdjf%ltotvel, sdjf%lzint, Kmm )
370      ELSE IF( LEN(TRIM(sdjf%wgtname)) > 0 ) THEN   ! On-the-fly interpolation
371         CALL wgt_list( sdjf, iw )
372         CALL fld_interp( sdjf%num, sdjf%clvar, iw, ipk, dta_alias(:,:,:), sdjf%nrec(1,iaa), sdjf%lsmname )
373         CALL lbc_lnk( 'fldread', dta_alias(:,:,:), sdjf%cltype, sdjf%zsgn, kfillmode = jpfillcopy )
374      ELSE                                          ! default case
375         ! C1D case: If product of spatial dimensions == ipk, then x,y are of
376         ! size 1 (point/mooring data): this must be read onto the central grid point
377         idvar  = iom_varid( sdjf%num, sdjf%clvar )
378         idmspc = iom_file ( sdjf%num )%ndims( idvar )
379         IF( iom_file( sdjf%num )%luld( idvar ) )   idmspc = idmspc - 1   ! id of the last spatial dimension
380         lmoor  = (  idmspc == 0 .OR. PRODUCT( iom_file( sdjf%num )%dimsz( 1:MAX(idmspc,1) ,idvar ) ) == ipk  )   
381         !
382         IF( lk_c1d .AND. lmoor ) THEN
383            CALL iom_get( sdjf%num, jpdom_unknown, sdjf%clvar, dta_alias(2,2,:), sdjf%nrec(1,iaa) )   ! jpdom_unknown -> no lbc_lnk
384            CALL lbc_lnk( 'fldread', dta_alias(:,:,:), 'T', 1., kfillmode = jpfillcopy )
385         ELSE
386            CALL iom_get( sdjf%num,  jpdom_global, sdjf%clvar, dta_alias(:,:,:), sdjf%nrec(1,iaa),   &
387               &          sdjf%cltype, sdjf%zsgn, kfill = jpfillcopy )
388         ENDIF
389      ENDIF
390      !
391      sdjf%rotn(iaa) = .false.   ! vector not yet rotated
392      !
393   END SUBROUTINE fld_get
394
395   
396   SUBROUTINE fld_map( knum, cdvar, pdta, krec, kmap, kgrd, kbdy, ldtotvel, ldzint, Kmm )
397      !!---------------------------------------------------------------------
398      !!                    ***  ROUTINE fld_map  ***
399      !!
400      !! ** Purpose :   read global data from file and map onto local data
401      !!                using a general mapping (for open boundaries)
402      !!----------------------------------------------------------------------
403      INTEGER                   , INTENT(in   ) ::   knum         ! stream number
404      CHARACTER(LEN=*)          , INTENT(in   ) ::   cdvar        ! variable name
405      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pdta         ! bdy output field on model grid
406      INTEGER                   , INTENT(in   ) ::   krec         ! record number to read (ie time slice)
407      INTEGER , DIMENSION(:)    , INTENT(in   ) ::   kmap         ! global-to-local bdy mapping indices
408      ! optional variables used for vertical interpolation:
409      INTEGER, OPTIONAL         , INTENT(in   ) ::   kgrd         ! grid type (t, u, v)
410      INTEGER, OPTIONAL         , INTENT(in   ) ::   kbdy         ! bdy number
411      LOGICAL, OPTIONAL         , INTENT(in   ) ::   ldtotvel     ! true if total ( = barotrop + barocline) velocity
412      LOGICAL, OPTIONAL         , INTENT(in   ) ::   ldzint       ! true if 3D variable requires a vertical interpolation
413      INTEGER, OPTIONAL         , INTENT(in   ) ::   Kmm          ! ocean time level index
414      !!
415      INTEGER                                   ::   ipi          ! length of boundary data on local process
416      INTEGER                                   ::   ipj          ! length of dummy dimension ( = 1 )
417      INTEGER                                   ::   ipk          ! number of vertical levels of pdta ( 2D: ipk=1 ; 3D: ipk=jpk )
418      INTEGER                                   ::   ipkb         ! number of vertical levels in boundary data file
419      INTEGER                                   ::   idvar        ! variable ID
420      INTEGER                                   ::   indims       ! number of dimensions of the variable
421      INTEGER, DIMENSION(4)                     ::   idimsz       ! size of variable dimensions
422      REAL(wp)                                  ::   zfv          ! fillvalue
423      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zz_read      ! work space for global boundary data
424      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zdta_read    ! work space local data requiring vertical interpolation
425      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zdta_read_z  ! work space local data requiring vertical interpolation
426      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zdta_read_dz ! work space local data requiring vertical interpolation
427      CHARACTER(LEN=1),DIMENSION(3)             ::   cltype
428      LOGICAL                                   ::   lluld        ! is the variable using the unlimited dimension
429      LOGICAL                                   ::   llzint       ! local value of ldzint
430      !!---------------------------------------------------------------------
431      !
432      cltype = (/'t','u','v'/)
433      !
434      ipi = SIZE( pdta, 1 )
435      ipj = SIZE( pdta, 2 )   ! must be equal to 1
436      ipk = SIZE( pdta, 3 )
437      !
438      llzint = .FALSE.
439      IF( PRESENT(ldzint) )   llzint = ldzint
440      !
441      idvar = iom_varid( knum, cdvar, kndims = indims, kdimsz = idimsz, lduld = lluld  )
442      IF( indims == 4 .OR. ( indims == 3 .AND. .NOT. lluld ) ) THEN   ;   ipkb = idimsz(3)   ! xy(zl)t or xy(zl)
443      ELSE                                                            ;   ipkb = 1           ! xy or xyt
444      ENDIF
445      !
446      ALLOCATE( zz_read( idimsz(1), idimsz(2), ipkb ) )  ! ++++++++ !!! this can be very big...         
447      !
448      IF( ipk == 1 ) THEN
449
450         IF( ipkb /= 1 ) CALL ctl_stop( 'fld_map : we must have ipkb = 1 to read surface data' )
451         CALL iom_get ( knum, jpdom_unknown, cdvar, zz_read(:,:,1), krec )   ! call iom_get with a 2D file
452         CALL fld_map_core( zz_read, kmap, pdta )
453
454      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
455      ! Do we include something here to adjust barotropic velocities !
456      ! in case of a depth difference between bdy files and          !
457      ! bathymetry in the case ln_totvel = .false. and ipkb>0?       !
458      ! [as the enveloping and parital cells could change H]         !
459      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
460
461      ELSE
462         !
463         CALL iom_get ( knum, jpdom_unknown, cdvar, zz_read(:,:,:), krec )   ! call iom_get with a 3D file
464         !
465         IF( ipkb /= ipk .OR. llzint ) THEN   ! boundary data not on model vertical grid : vertical interpolation
466            !
467            IF( ipk == jpk .AND. iom_varid(knum,'gdep'//cltype(kgrd)) /= -1 .AND. iom_varid(knum,'e3'//cltype(kgrd)) /= -1 ) THEN
468               
469               ALLOCATE( zdta_read(ipi,ipj,ipkb), zdta_read_z(ipi,ipj,ipkb), zdta_read_dz(ipi,ipj,ipkb) )
470               
471               CALL fld_map_core( zz_read, kmap, zdta_read )
472               CALL iom_get ( knum, jpdom_unknown, 'gdep'//cltype(kgrd), zz_read )   ! read only once? Potential temporal evolution?
473               CALL fld_map_core( zz_read, kmap, zdta_read_z )
474               CALL iom_get ( knum, jpdom_unknown,   'e3'//cltype(kgrd), zz_read )   ! read only once? Potential temporal evolution?
475               CALL fld_map_core( zz_read, kmap, zdta_read_dz )
476               
477               CALL iom_getatt(knum, '_FillValue', zfv, cdvar=cdvar )
478               CALL fld_bdy_interp(zdta_read, zdta_read_z, zdta_read_dz, pdta, kgrd, kbdy, zfv, ldtotvel, Kmm)
479               DEALLOCATE( zdta_read, zdta_read_z, zdta_read_dz )
480               
481            ELSE
482               IF( ipk /= jpk ) CALL ctl_stop( 'fld_map : this should be an impossible case...' )
483               WRITE(ctmp1,*) 'fld_map : vertical interpolation for bdy variable '//TRIM(cdvar)//' requires ' 
484               IF( iom_varid(knum, 'gdep'//cltype(kgrd)) == -1 ) CALL ctl_stop( ctmp1//'gdep'//cltype(kgrd)//' variable' )
485               IF( iom_varid(knum,   'e3'//cltype(kgrd)) == -1 ) CALL ctl_stop( ctmp1//  'e3'//cltype(kgrd)//' variable' )
486
487            ENDIF
488            !
489         ELSE                            ! bdy data assumed to be the same levels as bdy variables
490            !
491            CALL fld_map_core( zz_read, kmap, pdta )
492            !
493         ENDIF   ! ipkb /= ipk
494      ENDIF   ! ipk == 1
495     
496      DEALLOCATE( zz_read )
497
498   END SUBROUTINE fld_map
499
500     
501   SUBROUTINE fld_map_core( pdta_read, kmap, pdta_bdy )
502      !!---------------------------------------------------------------------
503      !!                    ***  ROUTINE fld_map_core  ***
504      !!
505      !! ** Purpose :  inner core of fld_map
506      !!----------------------------------------------------------------------
507      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read    ! global boundary data
508      INTEGER,  DIMENSION(:    ), INTENT(in   ) ::   kmap         ! global-to-local bdy mapping indices
509      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pdta_bdy     ! bdy output field on model grid
510      !!
511      INTEGER,  DIMENSION(3) ::   idim_read,  idim_bdy            ! arrays dimensions
512      INTEGER                ::   ji, jj, jk, jb                  ! loop counters
513      INTEGER                ::   im1
514      !!---------------------------------------------------------------------
515      !
516      idim_read = SHAPE( pdta_read )
517      idim_bdy  = SHAPE( pdta_bdy  )
518      !
519      ! in all cases: idim_bdy(2) == 1 .AND. idim_read(1) * idim_read(2) == idim_bdy(1)
520      ! structured BDY with rimwidth > 1                     : idim_read(2) == rimwidth /= 1
521      ! structured BDY with rimwidth == 1 or unstructured BDY: idim_read(2) == 1
522      !
523      IF( idim_read(2) > 1 ) THEN    ! structured BDY with rimwidth > 1 
524         DO jk = 1, idim_bdy(3)
525            DO jb = 1, idim_bdy(1)
526               im1 = kmap(jb) - 1
527               jj = im1 / idim_read(1) + 1
528               ji = MOD( im1, idim_read(1) ) + 1
529               pdta_bdy(jb,1,jk) =  pdta_read(ji,jj,jk)
530            END DO
531         END DO
532      ELSE
533         DO jk = 1, idim_bdy(3)
534            DO jb = 1, idim_bdy(1)   ! horizontal remap of bdy data on the local bdy
535               pdta_bdy(jb,1,jk) = pdta_read(kmap(jb),1,jk)
536            END DO
537         END DO
538      ENDIF
539     
540   END SUBROUTINE fld_map_core
541   
542   SUBROUTINE fld_bdy_interp(pdta_read, pdta_read_z, pdta_read_dz, pdta, kgrd, kbdy, pfv, ldtotvel, Kmm )
543      !!---------------------------------------------------------------------
544      !!                    ***  ROUTINE fld_bdy_interp  ***
545      !!
546      !! ** Purpose :   on the fly vertical interpolation to allow the use of
547      !!                boundary data from non-native vertical grid
548      !!----------------------------------------------------------------------
549      USE bdy_oce, ONLY:  idx_bdy         ! indexing for map <-> ij transformation
550
551      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read       ! data read in bdy file
552      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read_z     ! depth of the data read in bdy file
553      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read_dz    ! thickness of the levels in bdy file
554      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pdta            ! output field on model grid (2 dimensional)
555      REAL(wp)                  , INTENT(in   ) ::   pfv             ! fillvalue of the data read in bdy file
556      LOGICAL                   , INTENT(in   ) ::   ldtotvel        ! true if toal ( = barotrop + barocline) velocity
557      INTEGER                   , INTENT(in   ) ::   kgrd            ! grid type (t, u, v)
558      INTEGER                   , INTENT(in   ) ::   kbdy            ! bdy number
559      INTEGER, OPTIONAL         , INTENT(in   ) ::   Kmm             ! ocean time level index
560      !!
561      INTEGER                  ::   ipi                 ! length of boundary data on local process
562      INTEGER                  ::   ipkb                ! number of vertical levels in boundary data file
563      INTEGER                  ::   ipkmax              ! number of vertical levels in boundary data file where no mask
564      INTEGER                  ::   jb, ji, jj, jk, jkb ! loop counters
565      REAL(wp)                 ::   zcoef, zi           !
566      REAL(wp)                 ::   ztrans, ztrans_new  ! transports
567      REAL(wp), DIMENSION(jpk) ::   zdepth, zdhalf      ! level and half-level depth
568      !!---------------------------------------------------------------------
569     
570      ipi  = SIZE( pdta, 1 )
571      ipkb = SIZE( pdta_read, 3 )
572     
573      DO jb = 1, ipi
574         ji = idx_bdy(kbdy)%nbi(jb,kgrd)
575         jj = idx_bdy(kbdy)%nbj(jb,kgrd)
576         !
577         ! --- max jk where input data /= FillValue --- !
578         ipkmax = 1
579         DO jkb = 2, ipkb
580            IF( pdta_read(jb,1,jkb) /= pfv )   ipkmax = MAX( ipkmax, jkb )
581         END DO
582         !
583         ! --- calculate depth at t,u,v points --- !
584         SELECT CASE( kgrd )                         
585         CASE(1)            ! depth of T points:
586            zdepth(:) = gdept(ji,jj,:,Kmm)
587         CASE(2)            ! depth of U points: we must not use gdept_n as we don't want to do a communication
588            !                 --> copy what is done for gdept_n in domvvl...
589            zdhalf(1) = 0.0_wp
590            zdepth(1) = 0.5_wp * e3uw(ji,jj,1,Kmm)
591            DO jk = 2, jpk                               ! vertical sum
592               !    zcoef = umask - wumask    ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
593               !                              ! 1 everywhere from mbkt to mikt + 1 or 1 (if no isf)
594               !                              ! 0.5 where jk = mikt     
595               !!gm ???????   BUG ?  gdept_n as well as gde3w_n  does not include the thickness of ISF ??
596               zcoef = ( umask(ji,jj,jk) - wumask(ji,jj,jk) )
597               zdhalf(jk) = zdhalf(jk-1) + e3u(ji,jj,jk-1,Kmm)
598               zdepth(jk) =          zcoef  * ( zdhalf(jk  ) + 0.5_wp * e3uw(ji,jj,jk,Kmm))  &
599                  &         + (1._wp-zcoef) * ( zdepth(jk-1) +          e3uw(ji,jj,jk,Kmm))
600            END DO
601         CASE(3)            ! depth of V points: we must not use gdept_n as we don't want to do a communication
602            !                 --> copy what is done for gdept_n in domvvl...
603            zdhalf(1) = 0.0_wp
604            zdepth(1) = 0.5_wp * e3vw(ji,jj,1,Kmm)
605            DO jk = 2, jpk                               ! vertical sum
606               !    zcoef = vmask - wvmask    ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
607               !                              ! 1 everywhere from mbkt to mikt + 1 or 1 (if no isf)
608               !                              ! 0.5 where jk = mikt     
609               !!gm ???????   BUG ?  gdept_n as well as gde3w_n  does not include the thickness of ISF ??
610               zcoef = ( vmask(ji,jj,jk) - wvmask(ji,jj,jk) )
611               zdhalf(jk) = zdhalf(jk-1) + e3v(ji,jj,jk-1,Kmm)
612               zdepth(jk) =          zcoef  * ( zdhalf(jk  ) + 0.5_wp * e3vw(ji,jj,jk,Kmm))  &
613                  &         + (1._wp-zcoef) * ( zdepth(jk-1) +          e3vw(ji,jj,jk,Kmm))
614            END DO
615         END SELECT
616         !         
617         ! --- interpolate bdy data on the model grid --- !
618         DO jk = 1, jpk
619            IF(     zdepth(jk) <= pdta_read_z(jb,1,1)      ) THEN   ! above the first level of external data
620               pdta(jb,1,jk) = pdta_read(jb,1,1)
621            ELSEIF( zdepth(jk) >  pdta_read_z(jb,1,ipkmax) ) THEN   ! below the last level of external data /= FillValue
622               pdta(jb,1,jk) = pdta_read(jb,1,ipkmax)
623            ELSE                                                    ! inbetween: vertical interpolation between jkb & jkb+1
624               DO jkb = 1, ipkmax-1
625                  IF( ( ( zdepth(jk) - pdta_read_z(jb,1,jkb) ) * ( zdepth(jk) - pdta_read_z(jb,1,jkb+1) ) ) <= 0._wp ) THEN ! linear interpolation between 2 levels
626                     zi = ( zdepth(jk) - pdta_read_z(jb,1,jkb) ) / ( pdta_read_z(jb,1,jkb+1) - pdta_read_z(jb,1,jkb) )
627                     pdta(jb,1,jk) = pdta_read(jb,1,jkb) + zi * ( pdta_read(jb,1,jkb+1) - pdta_read(jb,1,jkb) )
628                  ENDIF
629               END DO
630            ENDIF
631         END DO   ! jpk
632         !
633      END DO   ! ipi
634
635      ! --- mask data and adjust transport --- !
636      SELECT CASE( kgrd )                         
637
638      CASE(1)                                 ! mask data (probably unecessary)
639         DO jb = 1, ipi
640            ji = idx_bdy(kbdy)%nbi(jb,kgrd)
641            jj = idx_bdy(kbdy)%nbj(jb,kgrd)
642            DO jk = 1, jpk                     
643               pdta(jb,1,jk) = pdta(jb,1,jk) * tmask(ji,jj,jk)
644            END DO
645         END DO
646         
647      CASE(2)                                 ! adjust the U-transport term
648         DO jb = 1, ipi
649            ji = idx_bdy(kbdy)%nbi(jb,kgrd)
650            jj = idx_bdy(kbdy)%nbj(jb,kgrd)
651            ztrans = 0._wp
652            DO jkb = 1, ipkb                              ! calculate transport on input grid
653               IF( pdta_read(jb,1,jkb) /= pfv )   ztrans = ztrans + pdta_read(jb,1,jkb) * pdta_read_dz(jb,1,jkb)
654            ENDDO
655            ztrans_new = 0._wp
656            DO jk = 1, jpk                                ! calculate transport on model grid
657               ztrans_new = ztrans_new +      pdta(jb,1,jk ) * e3u(ji,jj,jk,Kmm ) * umask(ji,jj,jk)
658            ENDDO
659            DO jk = 1, jpk                                ! make transport correction
660               IF(ldtotvel) THEN ! bdy data are total velocity so adjust bt transport term to match input data
661                  pdta(jb,1,jk) = ( pdta(jb,1,jk) + ( ztrans - ztrans_new ) * r1_hu(ji,jj,Kmm) ) * umask(ji,jj,jk)
662               ELSE              ! we're just dealing with bc velocity so bt transport term should sum to zero
663                  pdta(jb,1,jk) =   pdta(jb,1,jk) + (  0._wp - ztrans_new ) * r1_hu(ji,jj,Kmm)   * umask(ji,jj,jk)
664               ENDIF
665            ENDDO
666         ENDDO
667
668      CASE(3)                                 ! adjust the V-transport term
669         DO jb = 1, ipi
670            ji = idx_bdy(kbdy)%nbi(jb,kgrd)
671            jj = idx_bdy(kbdy)%nbj(jb,kgrd)
672            ztrans = 0._wp
673            DO jkb = 1, ipkb                              ! calculate transport on input grid
674               IF( pdta_read(jb,1,jkb) /= pfv )   ztrans = ztrans + pdta_read(jb,1,jkb) * pdta_read_dz(jb,1,jkb)
675            ENDDO
676            ztrans_new = 0._wp
677            DO jk = 1, jpk                                ! calculate transport on model grid
678               ztrans_new = ztrans_new +      pdta(jb,1,jk ) * e3v(ji,jj,jk,Kmm ) * vmask(ji,jj,jk)
679            ENDDO
680            DO jk = 1, jpk                                ! make transport correction
681               IF(ldtotvel) THEN ! bdy data are total velocity so adjust bt transport term to match input data
682                  pdta(jb,1,jk) = ( pdta(jb,1,jk) + ( ztrans - ztrans_new ) * r1_hv(ji,jj,Kmm) ) * vmask(ji,jj,jk)
683               ELSE              ! we're just dealing with bc velocity so bt transport term should sum to zero
684                  pdta(jb,1,jk) =   pdta(jb,1,jk) + (  0._wp - ztrans_new ) * r1_hv(ji,jj,Kmm)   * vmask(ji,jj,jk)
685               ENDIF
686            ENDDO
687         ENDDO
688      END SELECT
689     
690   END SUBROUTINE fld_bdy_interp
691
692
693   SUBROUTINE fld_rot( kt, sd )
694      !!---------------------------------------------------------------------
695      !!                    ***  ROUTINE fld_rot  ***
696      !!
697      !! ** Purpose :   Vector fields may need to be rotated onto the local grid direction
698      !!----------------------------------------------------------------------
699      INTEGER                , INTENT(in   ) ::   kt   ! ocean time step
700      TYPE(FLD), DIMENSION(:), INTENT(inout) ::   sd   ! input field related variables
701      !
702      INTEGER ::   ju, jv, jk, jn  ! loop indices
703      INTEGER ::   imf             ! size of the structure sd
704      INTEGER ::   ill             ! character length
705      INTEGER ::   iv              ! indice of V component
706      CHARACTER (LEN=100)          ::   clcomp       ! dummy weight name
707      REAL(wp), DIMENSION(jpi,jpj) ::   utmp, vtmp   ! temporary arrays for vector rotation
708      REAL(wp), DIMENSION(:,:,:), POINTER ::   dta_u, dta_v    ! short cut
709      !!---------------------------------------------------------------------
710      !
711      !! (sga: following code should be modified so that pairs arent searched for each time
712      !
713      imf = SIZE( sd )
714      DO ju = 1, imf
715         IF( TRIM(sd(ju)%clrootname) == 'NOT USED' )   CYCLE
716         ill = LEN_TRIM( sd(ju)%vcomp )
717         DO jn = 2-COUNT((/sd(ju)%ln_tint/)), 2
718            IF( ill > 0 .AND. .NOT. sd(ju)%rotn(jn) ) THEN   ! find vector rotations required             
719               IF( sd(ju)%vcomp(1:1) == 'U' ) THEN      ! east-west component has symbolic name starting with 'U'
720                  ! look for the north-south component which has same symbolic name but with 'U' replaced with 'V'
721                  clcomp = 'V' // sd(ju)%vcomp(2:ill)   ! works even if ill == 1
722                  iv = -1
723                  DO jv = 1, imf
724                     IF( TRIM(sd(jv)%clrootname) == 'NOT USED' )   CYCLE
725                     IF( TRIM(sd(jv)%vcomp) == TRIM(clcomp) )   iv = jv
726                  END DO
727                  IF( iv > 0 ) THEN   ! fields ju and iv are two components which need to be rotated together
728                     IF( sd(ju)%ln_tint ) THEN   ;   dta_u => sd(ju)%fdta(:,:,:,jn)   ;   dta_v => sd(iv)%fdta(:,:,:,jn) 
729                     ELSE                        ;   dta_u => sd(ju)%fnow(:,:,:   )   ;   dta_v => sd(iv)%fnow(:,:,:   )
730                     ENDIF
731                     DO jk = 1, SIZE( sd(ju)%fnow, 3 )
732                        CALL rot_rep( dta_u(:,:,jk), dta_v(:,:,jk), 'T', 'en->i', utmp(:,:) )
733                        CALL rot_rep( dta_u(:,:,jk), dta_v(:,:,jk), 'T', 'en->j', vtmp(:,:) )
734                        dta_u(:,:,jk) = utmp(:,:)   ;   dta_v(:,:,jk) = vtmp(:,:)
735                     END DO
736                     sd(ju)%rotn(jn) = .TRUE.               ! vector was rotated
737                     IF( lwp .AND. kt == nit000 )   WRITE(numout,*)   &
738                        &   'fld_read: vector pair ('//TRIM(sd(ju)%clvar)//', '//TRIM(sd(iv)%clvar)//') rotated on to model grid'
739                  ENDIF
740               ENDIF
741            ENDIF
742         END DO
743       END DO
744      !
745   END SUBROUTINE fld_rot
746
747
748   SUBROUTINE fld_def( sdjf, ldprev, ldnext )
749      !!---------------------------------------------------------------------
750      !!                    ***  ROUTINE fld_def  ***
751      !!
752      !! ** Purpose :   define the record(s) of the file and its name
753      !!----------------------------------------------------------------------
754      TYPE(FLD)        , INTENT(inout) ::   sdjf       ! input field related variables
755      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldprev     !
756      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldnext     !
757      !
758      INTEGER  :: jt
759      INTEGER  :: idaysec               ! number of seconds in 1 day = NINT(rday)
760      INTEGER  :: iyr, imt, idy, isecwk
761      INTEGER  :: indexyr, indexmt
762      INTEGER  :: ireclast
763      INTEGER  :: ishift, istart
764      INTEGER, DIMENSION(2)  :: isave
765      REAL(wp) :: zfreqs
766      LOGICAL  :: llprev, llnext, llstop
767      LOGICAL  :: llprevmt, llprevyr
768      LOGICAL  :: llnextmt, llnextyr
769      !!----------------------------------------------------------------------
770      idaysec = NINT(rday)
771      !
772      IF( PRESENT(ldprev) ) THEN   ;   llprev = ldprev
773      ELSE                         ;   llprev = .FALSE.
774      ENDIF
775      IF( PRESENT(ldnext) ) THEN   ;   llnext = ldnext
776      ELSE                         ;   llnext = .FALSE.
777      ENDIF
778
779      ! current file parameters
780      IF( sdjf%clftyp(1:4) == 'week' ) THEN         ! find the day of the beginning of the current week
781         isecwk = ksec_week( sdjf%clftyp(6:8) )     ! seconds between the beginning of the week and half of current time step
782         llprevmt = isecwk > nsec_month             ! longer time since beginning of the current week than the current month
783         llprevyr = llprevmt .AND. nmonth == 1
784         iyr = nyear  - COUNT((/llprevyr/))
785         imt = nmonth - COUNT((/llprevmt/)) + 12 * COUNT((/llprevyr/))
786         idy = nday + nmonth_len(nmonth-1) * COUNT((/llprevmt/)) - isecwk / idaysec
787         isecwk = nsec_year - isecwk                ! seconds between 00h jan 1st of current year and current week beginning
788      ELSE
789         iyr = nyear
790         imt = nmonth
791         idy = nday
792         isecwk  = 0
793      ENDIF
794
795      ! previous file parameters
796      IF( llprev ) THEN
797         IF( sdjf%clftyp(1:4) == 'week'    ) THEN   ! find the day of the beginning of previous week
798            isecwk = isecwk + 7 * idaysec           ! seconds between the beginning of previous week and half of the time step
799            llprevmt = isecwk > nsec_month          ! longer time since beginning of the previous week than the current month
800            llprevyr = llprevmt .AND. nmonth == 1
801            iyr = nyear  - COUNT((/llprevyr/))
802            imt = nmonth - COUNT((/llprevmt/)) + 12 * COUNT((/llprevyr/))
803            idy = nday + nmonth_len(nmonth-1) * COUNT((/llprevmt/)) - isecwk / idaysec
804            isecwk = nsec_year - isecwk             ! seconds between 00h jan 1st of current year and previous week beginning
805         ELSE
806            idy = nday   - COUNT((/ sdjf%clftyp == 'daily'                 /))
807            imt = nmonth - COUNT((/ sdjf%clftyp == 'monthly' .OR. idy == 0 /))
808            iyr = nyear  - COUNT((/ sdjf%clftyp == 'yearly'  .OR. imt == 0 /))
809            IF( idy == 0 ) idy = nmonth_len(imt)
810            IF( imt == 0 ) imt = 12
811            isecwk = 0
812         ENDIF
813      ENDIF
814
815      ! next file parameters
816      IF( llnext ) THEN
817         IF( sdjf%clftyp(1:4) == 'week'    ) THEN   ! find the day of the beginning of next week
818            isecwk = 7 * idaysec - isecwk           ! seconds between half of the time step and the beginning of next week
819            llnextmt = isecwk > ( nmonth_len(nmonth)*idaysec - nsec_month )   ! larger than the seconds to the end of the month
820            llnextyr = llnextmt .AND. nmonth == 12
821            iyr = nyear  + COUNT((/llnextyr/))
822            imt = nmonth + COUNT((/llnextmt/)) - 12 * COUNT((/llnextyr/))
823            idy = nday - nmonth_len(nmonth) * COUNT((/llnextmt/)) + isecwk / idaysec + 1
824            isecwk = nsec_year + isecwk             ! seconds between 00h jan 1st of current year and next week beginning
825         ELSE
826            idy = nday   + COUNT((/ sdjf%clftyp == 'daily'                                 /))
827            imt = nmonth + COUNT((/ sdjf%clftyp == 'monthly' .OR. idy > nmonth_len(nmonth) /))
828            iyr = nyear  + COUNT((/ sdjf%clftyp == 'yearly'  .OR. imt == 13                /))
829            IF( idy > nmonth_len(nmonth) )   idy = 1
830            IF( imt == 13                )   imt = 1
831            isecwk = 0
832         ENDIF
833      ENDIF
834      !
835      ! find the last record to be read -> update sdjf%nreclast
836      indexyr = iyr - nyear + 1                 ! which  year are we looking for? previous(0), current(1) or next(2)?
837      indexmt = imt + 12 * ( indexyr - 1 )      ! which month are we looking for (relatively to current year)?
838      !
839      ! Last record to be read in the current file
840      ! Predefine the number of record in the file according of its type.
841      ! We could compare this number with the number of records in the file and make a stop if the 2 numbers do not match...
842      ! However this would be much less fexible (e.g. for tests) and will force to rewite input files according to nleapy...
843      IF    ( NINT(sdjf%freqh) == -12 ) THEN            ;   ireclast = 1    ! yearly mean: consider only 1 record
844      ELSEIF( NINT(sdjf%freqh) ==  -1 ) THEN                                ! monthly mean:
845         IF(     sdjf%clftyp      == 'monthly' ) THEN   ;   ireclast = 1    !  consider that the file has  1 record
846         ELSE                                           ;   ireclast = 12   !  consider that the file has 12 record
847         ENDIF
848      ELSE                                                                  ! higher frequency mean (in hours)
849         IF(     sdjf%clftyp      == 'monthly' ) THEN   ;   ireclast = NINT( 24. * REAL(nmonth_len(indexmt), wp) / sdjf%freqh )
850         ELSEIF( sdjf%clftyp(1:4) == 'week'    ) THEN   ;   ireclast = NINT( 24. * 7.                            / sdjf%freqh )
851         ELSEIF( sdjf%clftyp      == 'daily'   ) THEN   ;   ireclast = NINT( 24.                                 / sdjf%freqh )
852         ELSE                                           ;   ireclast = NINT( 24. * REAL( nyear_len(indexyr), wp) / sdjf%freqh )
853         ENDIF
854      ENDIF
855
856      sdjf%nreclast = ireclast
857      ! Allocate arrays for beginning/middle/end of each record (seconds since Jan. 1st 00h of nit000 year)
858      IF( ALLOCATED(sdjf%nrecsec) )   DEALLOCATE( sdjf%nrecsec )
859      ALLOCATE( sdjf%nrecsec( 0:ireclast ) )
860      !
861      IF    ( NINT(sdjf%freqh) == -12 ) THEN                                     ! yearly mean and yearly file
862         SELECT CASE( indexyr )
863         CASE(0)   ;   sdjf%nrecsec(0) = nsec1jan000 - nyear_len( 0 ) * idaysec
864         CASE(1)   ;   sdjf%nrecsec(0) = nsec1jan000
865         CASE(2)   ;   sdjf%nrecsec(0) = nsec1jan000 + nyear_len( 1 ) * idaysec
866         ENDSELECT
867         sdjf%nrecsec(1) = sdjf%nrecsec(0) + nyear_len( indexyr ) * idaysec
868      ELSEIF( NINT(sdjf%freqh) ==  -1 ) THEN                                     ! monthly mean:
869         IF(     sdjf%clftyp      == 'monthly' ) THEN                            !    monthly file
870            sdjf%nrecsec(0   ) = nsec1jan000 + nmonth_beg(indexmt  )
871            sdjf%nrecsec(1   ) = nsec1jan000 + nmonth_beg(indexmt+1)
872         ELSE                                                                    !    yearly  file
873            ishift = 12 * ( indexyr - 1 )
874            sdjf%nrecsec(0:12) = nsec1jan000 + nmonth_beg(1+ishift:13+ishift)
875         ENDIF
876      ELSE                                                                       ! higher frequency mean (in hours)
877         IF(     sdjf%clftyp      == 'monthly' ) THEN   ;   istart = nsec1jan000 + nmonth_beg(indexmt)
878         ELSEIF( sdjf%clftyp(1:4) == 'week'    ) THEN   ;   istart = nsec1jan000 + isecwk
879         ELSEIF( sdjf%clftyp      == 'daily'   ) THEN   ;   istart = nsec1jan000 + nmonth_beg(indexmt) + ( idy - 1 ) * idaysec
880         ELSEIF( indexyr          == 0         ) THEN   ;   istart = nsec1jan000 - nyear_len( 0 ) * idaysec
881         ELSEIF( indexyr          == 2         ) THEN   ;   istart = nsec1jan000 + nyear_len( 1 ) * idaysec
882         ELSE                                           ;   istart = nsec1jan000
883         ENDIF
884         zfreqs = sdjf%freqh * rhhmm * rmmss
885         DO jt = 0, sdjf%nreclast
886            sdjf%nrecsec(jt) = istart + NINT( zfreqs * REAL(jt,wp) )
887         END DO
888      ENDIF
889      !
890      IF( sdjf%ln_tint ) THEN   ! record time defined in the middle of the record, computed using an implementation
891                                ! of the rounded average that is valid over the full integer range
892         sdjf%nrecsec(1:sdjf%nreclast) = sdjf%nrecsec(0:sdjf%nreclast-1) / 2 + sdjf%nrecsec(1:sdjf%nreclast) / 2 + &
893            & MAX( MOD( sdjf%nrecsec(0:sdjf%nreclast-1), 2 ), MOD( sdjf%nrecsec(1:sdjf%nreclast), 2 ) )
894      END IF
895      !
896      sdjf%clname = fld_filename( sdjf, idy, imt, iyr )
897      !
898   END SUBROUTINE fld_def
899
900   
901   SUBROUTINE fld_clopn( sdjf )
902      !!---------------------------------------------------------------------
903      !!                    ***  ROUTINE fld_clopn  ***
904      !!
905      !! ** Purpose :   close/open the files
906      !!----------------------------------------------------------------------
907      TYPE(FLD)        , INTENT(inout) ::   sdjf       ! input field related variables
908      !
909      INTEGER, DIMENSION(2)  :: isave
910      LOGICAL  :: llprev, llnext, llstop
911      !!----------------------------------------------------------------------
912      !
913      llprev = sdjf%nrecsec(sdjf%nreclast) < nsec000_1jan000   ! file ends before the beginning of the job -> file may not exist
914      llnext = sdjf%nrecsec(       0     ) > nsecend_1jan000   ! file begins after the end of the job -> file may not exist
915
916      llstop = sdjf%ln_clim .OR. .NOT. ( llprev .OR. llnext )
917
918      IF( sdjf%num <= 0 .OR. .NOT. sdjf%ln_clim  ) THEN
919         IF( sdjf%num > 0 )   CALL iom_close( sdjf%num )   ! close file if already open
920         CALL iom_open( sdjf%clname, sdjf%num, ldstop = llstop, ldiof = LEN_TRIM(sdjf%wgtname) > 0 )
921      ENDIF
922      !
923      IF( sdjf%num <= 0 .AND. .NOT. llstop ) THEN   ! file not found but we do accept this...
924         !
925         IF( llprev ) THEN   ! previous file does not exist : go back to current and accept to read only the first record
926            CALL ctl_warn('previous file: '//TRIM(sdjf%clname)//' not found -> go back to current year/month/week/day file')
927            isave(1:2) = sdjf%nrecsec(sdjf%nreclast-1:sdjf%nreclast)   ! save previous file info
928            CALL fld_def( sdjf )   ! go back to current file
929            sdjf%nreclast = 1   ! force to use only the first record (do as if other were not existing...)
930            sdjf%nrecsec(0:1) = isave(1:2)
931         ENDIF
932         !
933         IF( llnext ) THEN   ! next     file does not exist : go back to current and accept to read only the last  record
934            CALL ctl_warn('next file: '//TRIM(sdjf%clname)//' not found -> go back to current year/month/week/day file')
935            isave(1:2) = sdjf%nrecsec(0:1)    ! save next file info
936            CALL fld_def( sdjf )   ! go back to current file
937            ! -> read last record but keep record info from the first record of next file
938            sdjf%nrecsec(sdjf%nreclast-1:sdjf%nreclast) = isave(1:2)
939            sdjf%nrecsec(0:sdjf%nreclast-2) = nflag
940         ENDIF
941         !
942         CALL iom_open( sdjf%clname, sdjf%num, ldiof = LEN_TRIM(sdjf%wgtname) > 0 )   
943         !
944      ENDIF
945      !
946   END SUBROUTINE fld_clopn
947
948
949   SUBROUTINE fld_fill( sdf, sdf_n, cdir, cdcaller, cdtitle, cdnam, knoprint )
950      !!---------------------------------------------------------------------
951      !!                    ***  ROUTINE fld_fill  ***
952      !!
953      !! ** Purpose :   fill the data structure (sdf) with the associated information
954      !!              read in namelist (sdf_n) and control print
955      !!----------------------------------------------------------------------
956      TYPE(FLD)  , DIMENSION(:)          , INTENT(inout) ::   sdf        ! structure of input fields (file informations, fields read)
957      TYPE(FLD_N), DIMENSION(:)          , INTENT(in   ) ::   sdf_n      ! array of namelist information structures
958      CHARACTER(len=*)                   , INTENT(in   ) ::   cdir       ! Root directory for location of flx files
959      CHARACTER(len=*)                   , INTENT(in   ) ::   cdcaller   ! name of the calling routine
960      CHARACTER(len=*)                   , INTENT(in   ) ::   cdtitle    ! description of the calling routine
961      CHARACTER(len=*)                   , INTENT(in   ) ::   cdnam      ! name of the namelist from which sdf_n comes
962      INTEGER                  , OPTIONAL, INTENT(in   ) ::   knoprint   ! no calling routine information printed
963      !
964      INTEGER  ::   jf   ! dummy indices
965      !!---------------------------------------------------------------------
966      !
967      DO jf = 1, SIZE(sdf)
968         sdf(jf)%clrootname = sdf_n(jf)%clname
969         IF( TRIM(sdf_n(jf)%clname) /= 'NOT USED' )   sdf(jf)%clrootname = TRIM( cdir )//sdf(jf)%clrootname
970         sdf(jf)%clname     = "not yet defined"
971         sdf(jf)%freqh      = sdf_n(jf)%freqh
972         sdf(jf)%clvar      = sdf_n(jf)%clvar
973         sdf(jf)%ln_tint    = sdf_n(jf)%ln_tint
974         sdf(jf)%ln_clim    = sdf_n(jf)%ln_clim
975         sdf(jf)%clftyp     = sdf_n(jf)%clftyp
976         sdf(jf)%cltype     = 'T'   ! by default don't do any call to lbc_lnk in iom_get
977         sdf(jf)%zsgn       = 1.    ! by default don't do change signe across the north fold
978         sdf(jf)%num        = -1
979         sdf(jf)%nbb        = 1  ! start with before data in 1
980         sdf(jf)%naa        = 2  ! start with after  data in 2
981         sdf(jf)%wgtname    = " "
982         IF( LEN( TRIM(sdf_n(jf)%wname) ) > 0 )   sdf(jf)%wgtname = TRIM( cdir )//sdf_n(jf)%wname
983         sdf(jf)%lsmname = " "
984         IF( LEN( TRIM(sdf_n(jf)%lname) ) > 0 )   sdf(jf)%lsmname = TRIM( cdir )//sdf_n(jf)%lname
985         sdf(jf)%vcomp      = sdf_n(jf)%vcomp
986         sdf(jf)%rotn(:)    = .TRUE.   ! pretend to be rotated -> won't try to rotate data before the first call to fld_get
987         IF( sdf(jf)%clftyp(1:4) == 'week' .AND. nn_leapy == 0  )   &
988            &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdf(jf)%clrootname)//') needs nn_leapy = 1')
989         IF( sdf(jf)%clftyp(1:4) == 'week' .AND. sdf(jf)%ln_clim )   &
990            &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdf(jf)%clrootname)//') needs ln_clim = .FALSE.')
991         sdf(jf)%nreclast   = -1 ! Set to non zero default value to avoid errors, is updated to meaningful value during fld_clopn
992         sdf(jf)%igrd       = 0
993         sdf(jf)%ibdy       = 0
994         sdf(jf)%imap       => NULL()
995         sdf(jf)%ltotvel    = .FALSE.
996         sdf(jf)%lzint      = .FALSE.
997      END DO
998      !
999      IF(lwp) THEN      ! control print
1000         WRITE(numout,*)
1001         IF( .NOT.PRESENT( knoprint) ) THEN
1002            WRITE(numout,*) TRIM( cdcaller )//' : '//TRIM( cdtitle )
1003            WRITE(numout,*) (/ ('~', jf = 1, LEN_TRIM( cdcaller ) ) /)
1004         ENDIF
1005         WRITE(numout,*) '   fld_fill : fill data structure with information from namelist '//TRIM( cdnam )
1006         WRITE(numout,*) '   ~~~~~~~~'
1007         WRITE(numout,*) '      list of files and frequency (>0: in hours ; <0 in months)'
1008         DO jf = 1, SIZE(sdf)
1009            WRITE(numout,*) '      root filename: '  , TRIM( sdf(jf)%clrootname ), '   variable name: ', TRIM( sdf(jf)%clvar )
1010            WRITE(numout,*) '         frequency: '      ,       sdf(jf)%freqh       ,   &
1011               &                  '   time interp: '    ,       sdf(jf)%ln_tint     ,   &
1012               &                  '   climatology: '    ,       sdf(jf)%ln_clim
1013            WRITE(numout,*) '         weights: '        , TRIM( sdf(jf)%wgtname    ),   &
1014               &                  '   pairing: '        , TRIM( sdf(jf)%vcomp      ),   &
1015               &                  '   data type: '      ,       sdf(jf)%clftyp      ,   &
1016               &                  '   land/sea mask:'   , TRIM( sdf(jf)%lsmname    )
1017            call flush(numout)
1018         END DO
1019      ENDIF
1020      !
1021   END SUBROUTINE fld_fill
1022
1023
1024   SUBROUTINE wgt_list( sd, kwgt )
1025      !!---------------------------------------------------------------------
1026      !!                    ***  ROUTINE wgt_list  ***
1027      !!
1028      !! ** Purpose :   search array of WGTs and find a weights file entry,
1029      !!                or return a new one adding it to the end if new entry.
1030      !!                the weights data is read in and restructured (fld_weight)
1031      !!----------------------------------------------------------------------
1032      TYPE( FLD ), INTENT(in   ) ::   sd        ! field with name of weights file
1033      INTEGER    , INTENT(  out) ::   kwgt      ! index of weights
1034      !
1035      INTEGER ::   kw, nestid   ! local integer
1036      !!----------------------------------------------------------------------
1037      !
1038      !! search down linked list
1039      !! weights filename is either present or we hit the end of the list
1040      !
1041      !! because agrif nest part of filenames are now added in iom_open
1042      !! to distinguish between weights files on the different grids, need to track
1043      !! nest number explicitly
1044      nestid = 0
1045#if defined key_agrif
1046      nestid = Agrif_Fixed()
1047#endif
1048      DO kw = 1, nxt_wgt-1
1049         IF( ref_wgts(kw)%wgtname == sd%wgtname .AND. &
1050             ref_wgts(kw)%nestid  == nestid) THEN
1051            kwgt = kw
1052            RETURN
1053         ENDIF
1054      END DO
1055      kwgt = nxt_wgt
1056      CALL fld_weight( sd )
1057      !
1058   END SUBROUTINE wgt_list
1059
1060
1061   SUBROUTINE wgt_print( )
1062      !!---------------------------------------------------------------------
1063      !!                    ***  ROUTINE wgt_print  ***
1064      !!
1065      !! ** Purpose :   print the list of known weights
1066      !!----------------------------------------------------------------------
1067      INTEGER ::   kw   !
1068      !!----------------------------------------------------------------------
1069      !
1070      DO kw = 1, nxt_wgt-1
1071         WRITE(numout,*) 'weight file:  ',TRIM(ref_wgts(kw)%wgtname)
1072         WRITE(numout,*) '      ddims:  ',ref_wgts(kw)%ddims(1),ref_wgts(kw)%ddims(2)
1073         WRITE(numout,*) '     numwgt:  ',ref_wgts(kw)%numwgt
1074         WRITE(numout,*) '     jpiwgt:  ',ref_wgts(kw)%jpiwgt
1075         WRITE(numout,*) '     jpjwgt:  ',ref_wgts(kw)%jpjwgt
1076         WRITE(numout,*) '    botleft:  ',ref_wgts(kw)%botleft
1077         WRITE(numout,*) '   topright:  ',ref_wgts(kw)%topright
1078         IF( ref_wgts(kw)%cyclic ) THEN
1079            WRITE(numout,*) '       cyclical'
1080            IF( ref_wgts(kw)%overlap > 0 ) WRITE(numout,*) '              with overlap of ', ref_wgts(kw)%overlap
1081         ELSE
1082            WRITE(numout,*) '       not cyclical'
1083         ENDIF
1084         IF( ASSOCIATED(ref_wgts(kw)%data_wgt) )  WRITE(numout,*) '       allocated'
1085      END DO
1086      !
1087   END SUBROUTINE wgt_print
1088
1089
1090   SUBROUTINE fld_weight( sd )
1091      !!---------------------------------------------------------------------
1092      !!                    ***  ROUTINE fld_weight  ***
1093      !!
1094      !! ** Purpose :   create a new WGT structure and fill in data from file,
1095      !!              restructuring as required
1096      !!----------------------------------------------------------------------
1097      TYPE( FLD ), INTENT(in) ::   sd   ! field with name of weights file
1098      !!
1099      INTEGER ::   ji,jj,jn   ! dummy loop indices
1100      INTEGER ::   inum       ! local logical unit
1101      INTEGER ::   id         ! local variable id
1102      INTEGER ::   ipk        ! local vertical dimension
1103      INTEGER ::   zwrap      ! local integer
1104      LOGICAL ::   cyclical   !
1105      CHARACTER (len=5) ::   clname   !
1106      INTEGER , DIMENSION(4) ::   ddims
1107      INTEGER                ::   isrc
1108      REAL(wp), DIMENSION(jpi,jpj) ::   data_tmp
1109      !!----------------------------------------------------------------------
1110      !
1111      IF( nxt_wgt > tot_wgts ) THEN
1112        CALL ctl_stop("fld_weight: weights array size exceeded, increase tot_wgts")
1113      ENDIF
1114      !
1115      !! new weights file entry, add in extra information
1116      !! a weights file represents a 2D grid of a certain shape, so we assume that the current
1117      !! input data file is representative of all other files to be opened and processed with the
1118      !! current weights file
1119
1120      !! get data grid dimensions
1121      id = iom_varid( sd%num, sd%clvar, ddims )
1122
1123      !! now open the weights file
1124      CALL iom_open ( sd%wgtname, inum )   ! interpolation weights
1125      IF( inum > 0 ) THEN
1126
1127         !! determine whether we have an east-west cyclic grid
1128         !! from global attribute called "ew_wrap" in the weights file
1129         !! note that if not found, iom_getatt returns -999 and cyclic with no overlap is assumed
1130         !! since this is the most common forcing configuration
1131
1132         CALL iom_getatt(inum, 'ew_wrap', zwrap)
1133         IF( zwrap >= 0 ) THEN
1134            cyclical = .TRUE.
1135         ELSE IF( zwrap == -999 ) THEN
1136            cyclical = .TRUE.
1137            zwrap = 0
1138         ELSE
1139            cyclical = .FALSE.
1140         ENDIF
1141
1142         ref_wgts(nxt_wgt)%ddims(1) = ddims(1)
1143         ref_wgts(nxt_wgt)%ddims(2) = ddims(2)
1144         ref_wgts(nxt_wgt)%wgtname = sd%wgtname
1145         ref_wgts(nxt_wgt)%overlap = zwrap
1146         ref_wgts(nxt_wgt)%cyclic = cyclical
1147         ref_wgts(nxt_wgt)%nestid = 0
1148#if defined key_agrif
1149         ref_wgts(nxt_wgt)%nestid = Agrif_Fixed()
1150#endif
1151         !! weights file is stored as a set of weights (wgt01->wgt04 or wgt01->wgt16)
1152         !! for each weight wgtNN there is an integer array srcNN which gives the point in
1153         !! the input data grid which is to be multiplied by the weight
1154         !! they are both arrays on the model grid so the result of the multiplication is
1155         !! added into an output array on the model grid as a running sum
1156
1157         !! two possible cases: bilinear (4 weights) or bicubic (16 weights)
1158         id = iom_varid(inum, 'src05', ldstop=.FALSE.)
1159         IF( id <= 0 ) THEN   ;   ref_wgts(nxt_wgt)%numwgt = 4
1160         ELSE                 ;   ref_wgts(nxt_wgt)%numwgt = 16
1161         ENDIF
1162
1163         ALLOCATE( ref_wgts(nxt_wgt)%data_jpi(Nis0:Nie0,Njs0:Nje0,4) )
1164         ALLOCATE( ref_wgts(nxt_wgt)%data_jpj(Nis0:Nie0,Njs0:Nje0,4) )
1165         ALLOCATE( ref_wgts(nxt_wgt)%data_wgt(Nis0:Nie0,Njs0:Nje0,ref_wgts(nxt_wgt)%numwgt) )
1166
1167         DO jn = 1,4
1168            WRITE(clname,'(a3,i2.2)') 'src',jn
1169            CALL iom_get ( inum, jpdom_global, clname, data_tmp(:,:), cd_type = 'Z' )   !  no call to lbc_lnk
1170            DO_2D_00_00
1171               isrc = NINT(data_tmp(ji,jj)) - 1
1172               ref_wgts(nxt_wgt)%data_jpi(ji,jj,jn) = 1 + MOD(isrc,  ref_wgts(nxt_wgt)%ddims(1))
1173               ref_wgts(nxt_wgt)%data_jpj(ji,jj,jn) = 1 +     isrc / ref_wgts(nxt_wgt)%ddims(1)
1174            END_2D
1175         END DO
1176
1177         DO jn = 1, ref_wgts(nxt_wgt)%numwgt
1178            WRITE(clname,'(a3,i2.2)') 'wgt',jn
1179            CALL iom_get ( inum, jpdom_global, clname, data_tmp(:,:), cd_type = 'Z' )   !  no call to lbc_lnk
1180            DO_2D_00_00
1181               ref_wgts(nxt_wgt)%data_wgt(ji,jj,jn) = data_tmp(ji,jj)
1182            END_2D
1183         END DO
1184         CALL iom_close (inum)
1185 
1186         ! find min and max indices in grid
1187         ref_wgts(nxt_wgt)%botleft( 1) = MINVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1188         ref_wgts(nxt_wgt)%botleft( 2) = MINVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1189         ref_wgts(nxt_wgt)%topright(1) = MAXVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1190         ref_wgts(nxt_wgt)%topright(2) = MAXVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1191
1192         ! and therefore dimensions of the input box
1193         ref_wgts(nxt_wgt)%jpiwgt = ref_wgts(nxt_wgt)%topright(1) - ref_wgts(nxt_wgt)%botleft(1) + 1
1194         ref_wgts(nxt_wgt)%jpjwgt = ref_wgts(nxt_wgt)%topright(2) - ref_wgts(nxt_wgt)%botleft(2) + 1
1195
1196         ! shift indexing of source grid
1197         ref_wgts(nxt_wgt)%data_jpi(:,:,:) = ref_wgts(nxt_wgt)%data_jpi(:,:,:) - ref_wgts(nxt_wgt)%botleft(1) + 1
1198         ref_wgts(nxt_wgt)%data_jpj(:,:,:) = ref_wgts(nxt_wgt)%data_jpj(:,:,:) - ref_wgts(nxt_wgt)%botleft(2) + 1
1199
1200         ! create input grid, give it a halo to allow gradient calculations
1201         ! SA: +3 stencil is a patch to avoid out-of-bound computation in some configuration.
1202         ! a more robust solution will be given in next release
1203         ipk =  SIZE(sd%fnow, 3)
1204         ALLOCATE( ref_wgts(nxt_wgt)%fly_dta(ref_wgts(nxt_wgt)%jpiwgt+3, ref_wgts(nxt_wgt)%jpjwgt+3 ,ipk) )
1205         IF( ref_wgts(nxt_wgt)%cyclic ) ALLOCATE( ref_wgts(nxt_wgt)%col(1,ref_wgts(nxt_wgt)%jpjwgt+3,ipk) )
1206         !
1207         nxt_wgt = nxt_wgt + 1
1208         !
1209      ELSE
1210         CALL ctl_stop( '    fld_weight : unable to read the file ' )
1211      ENDIF
1212      !
1213   END SUBROUTINE fld_weight
1214
1215
1216   SUBROUTINE apply_seaoverland( clmaskfile, zfieldo, jpi1_lsm, jpi2_lsm, jpj1_lsm,   &
1217      &                          jpj2_lsm, itmpi, itmpj, itmpz, rec1_lsm, recn_lsm )
1218      !!---------------------------------------------------------------------
1219      !!                    ***  ROUTINE apply_seaoverland  ***
1220      !!
1221      !! ** Purpose :   avoid spurious fluxes in coastal or near-coastal areas
1222      !!                due to the wrong usage of "land" values from the coarse
1223      !!                atmospheric model when spatial interpolation is required
1224      !!      D. Delrosso INGV         
1225      !!----------------------------------------------------------------------
1226      INTEGER,                   INTENT(in   ) :: itmpi,itmpj,itmpz                    ! lengths
1227      INTEGER,                   INTENT(in   ) :: jpi1_lsm,jpi2_lsm,jpj1_lsm,jpj2_lsm  ! temporary indices
1228      INTEGER, DIMENSION(3),     INTENT(in   ) :: rec1_lsm,recn_lsm                    ! temporary arrays for start and length
1229      REAL(wp),DIMENSION (:,:,:),INTENT(inout) :: zfieldo                              ! input/output array for seaoverland application
1230      CHARACTER (len=100),       INTENT(in   ) :: clmaskfile                           ! land/sea mask file name
1231      !
1232      INTEGER :: inum,jni,jnj,jnz,jc   ! local indices
1233      REAL(wp),DIMENSION (:,:,:),ALLOCATABLE :: zslmec1             ! local array for land point detection
1234      REAL(wp),DIMENSION (:,:),  ALLOCATABLE :: zfieldn   ! array of forcing field with undeff for land points
1235      REAL(wp),DIMENSION (:,:),  ALLOCATABLE :: zfield    ! array of forcing field
1236      !!---------------------------------------------------------------------
1237      !
1238      ALLOCATE ( zslmec1(itmpi,itmpj,itmpz), zfieldn(itmpi,itmpj), zfield(itmpi,itmpj) )
1239      !
1240      ! Retrieve the land sea mask data
1241      CALL iom_open( clmaskfile, inum )
1242      SELECT CASE( SIZE(zfieldo(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),3) )
1243      CASE(1)
1244         CALL iom_get( inum, jpdom_unknown, 'LSM', zslmec1(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,1),   &
1245            &          1, kstart = rec1_lsm, kcount = recn_lsm)
1246      CASE DEFAULT
1247         CALL iom_get( inum, jpdom_unknown, 'LSM', zslmec1(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),   &
1248            &          1, kstart = rec1_lsm, kcount = recn_lsm)
1249      END SELECT
1250      CALL iom_close( inum )
1251      !
1252      DO jnz=1,rec1_lsm(3)             !! Loop over k dimension
1253         !
1254         DO jni = 1, itmpi                               !! copy the original field into a tmp array
1255            DO jnj = 1, itmpj                            !! substituting undeff over land points
1256               zfieldn(jni,jnj) = zfieldo(jni,jnj,jnz)
1257               IF( zslmec1(jni,jnj,jnz) == 1. )   zfieldn(jni,jnj) = undeff_lsm
1258            END DO
1259         END DO
1260         !
1261         CALL seaoverland( zfieldn, itmpi, itmpj, zfield )
1262         DO jc = 1, nn_lsm
1263            CALL seaoverland( zfield, itmpi, itmpj, zfield )
1264         END DO
1265         !
1266         !   Check for Undeff and substitute original values
1267         IF( ANY(zfield==undeff_lsm) ) THEN
1268            DO jni = 1, itmpi
1269               DO jnj = 1, itmpj
1270                  IF( zfield(jni,jnj)==undeff_lsm )   zfield(jni,jnj) = zfieldo(jni,jnj,jnz)
1271               END DO
1272            END DO
1273         ENDIF
1274         !
1275         zfieldo(:,:,jnz) = zfield(:,:)
1276         !
1277      END DO                           !! End Loop over k dimension
1278      !
1279      DEALLOCATE ( zslmec1, zfieldn, zfield )
1280      !
1281   END SUBROUTINE apply_seaoverland 
1282
1283
1284   SUBROUTINE seaoverland( zfieldn, ileni, ilenj, zfield )
1285      !!---------------------------------------------------------------------
1286      !!                    ***  ROUTINE seaoverland  ***
1287      !!
1288      !! ** Purpose :   create shifted matrices for seaoverland application 
1289      !!      D. Delrosso INGV
1290      !!----------------------------------------------------------------------
1291      INTEGER                      , INTENT(in   ) :: ileni,ilenj   ! lengths
1292      REAL, DIMENSION (ileni,ilenj), INTENT(in   ) :: zfieldn       ! array of forcing field with undeff for land points
1293      REAL, DIMENSION (ileni,ilenj), INTENT(  out) :: zfield        ! array of forcing field
1294      !
1295      REAL   , DIMENSION (ileni,ilenj)   :: zmat1, zmat2, zmat3, zmat4  ! local arrays
1296      REAL   , DIMENSION (ileni,ilenj)   :: zmat5, zmat6, zmat7, zmat8  !   -     -
1297      REAL   , DIMENSION (ileni,ilenj)   :: zlsm2d                      !   -     -
1298      REAL   , DIMENSION (ileni,ilenj,8) :: zlsm3d                      !   -     -
1299      LOGICAL, DIMENSION (ileni,ilenj,8) :: ll_msknan3d                 ! logical mask for undeff detection
1300      LOGICAL, DIMENSION (ileni,ilenj)   :: ll_msknan2d                 ! logical mask for undeff detection
1301      !!----------------------------------------------------------------------
1302      zmat8 = eoshift( zfieldn , SHIFT=-1 , BOUNDARY = (/zfieldn(:,1)/)     , DIM=2 )
1303      zmat1 = eoshift( zmat8   , SHIFT=-1 , BOUNDARY = (/zmat8(1,:)/)       , DIM=1 )
1304      zmat2 = eoshift( zfieldn , SHIFT=-1 , BOUNDARY = (/zfieldn(1,:)/)     , DIM=1 )
1305      zmat4 = eoshift( zfieldn , SHIFT= 1 , BOUNDARY = (/zfieldn(:,ilenj)/) , DIM=2 )
1306      zmat3 = eoshift( zmat4   , SHIFT=-1 , BOUNDARY = (/zmat4(1,:)/)       , DIM=1 )
1307      zmat5 = eoshift( zmat4   , SHIFT= 1 , BOUNDARY = (/zmat4(ileni,:)/)   , DIM=1 )
1308      zmat6 = eoshift( zfieldn , SHIFT= 1 , BOUNDARY = (/zfieldn(ileni,:)/) , DIM=1 )
1309      zmat7 = eoshift( zmat8   , SHIFT= 1 , BOUNDARY = (/zmat8(ileni,:)/)   , DIM=1 )
1310      !
1311      zlsm3d  = RESHAPE( (/ zmat1, zmat2, zmat3, zmat4, zmat5, zmat6, zmat7, zmat8 /), (/ ileni, ilenj, 8 /))
1312      ll_msknan3d = .NOT.( zlsm3d  == undeff_lsm )
1313      ll_msknan2d = .NOT.( zfieldn == undeff_lsm )  ! FALSE where is Undeff (land)
1314      zlsm2d = SUM( zlsm3d, 3 , ll_msknan3d ) / MAX( 1 , COUNT( ll_msknan3d , 3 ) )
1315      WHERE( COUNT( ll_msknan3d , 3 ) == 0._wp )   zlsm2d = undeff_lsm
1316      zfield = MERGE( zfieldn, zlsm2d, ll_msknan2d )
1317      !
1318   END SUBROUTINE seaoverland
1319
1320
1321   SUBROUTINE fld_interp( num, clvar, kw, kk, dta, nrec, lsmfile)     
1322      !!---------------------------------------------------------------------
1323      !!                    ***  ROUTINE fld_interp  ***
1324      !!
1325      !! ** Purpose :   apply weights to input gridded data to create data
1326      !!                on model grid
1327      !!----------------------------------------------------------------------
1328      INTEGER                   , INTENT(in   ) ::   num     ! stream number
1329      CHARACTER(LEN=*)          , INTENT(in   ) ::   clvar   ! variable name
1330      INTEGER                   , INTENT(in   ) ::   kw      ! weights number
1331      INTEGER                   , INTENT(in   ) ::   kk      ! vertical dimension of kk
1332      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   dta     ! output field on model grid
1333      INTEGER                   , INTENT(in   ) ::   nrec    ! record number to read (ie time slice)
1334      CHARACTER(LEN=*)          , INTENT(in   ) ::   lsmfile ! land sea mask file name
1335      !
1336      INTEGER, DIMENSION(3) ::   rec1, recn           ! temporary arrays for start and length
1337      INTEGER, DIMENSION(3) ::   rec1_lsm, recn_lsm   ! temporary arrays for start and length in case of seaoverland
1338      INTEGER ::   ii_lsm1,ii_lsm2,ij_lsm1,ij_lsm2    ! temporary indices
1339      INTEGER ::   ji, jj, jk, jn, jir, jjr           ! loop counters
1340      INTEGER ::   ipk
1341      INTEGER ::   ni, nj                             ! lengths
1342      INTEGER ::   jpimin,jpiwid                      ! temporary indices
1343      INTEGER ::   jpimin_lsm,jpiwid_lsm              ! temporary indices
1344      INTEGER ::   jpjmin,jpjwid                      ! temporary indices
1345      INTEGER ::   jpjmin_lsm,jpjwid_lsm              ! temporary indices
1346      INTEGER ::   jpi1,jpi2,jpj1,jpj2                ! temporary indices
1347      INTEGER ::   jpi1_lsm,jpi2_lsm,jpj1_lsm,jpj2_lsm   ! temporary indices
1348      INTEGER ::   itmpi,itmpj,itmpz                     ! lengths
1349      REAL(wp),DIMENSION(:,:,:), ALLOCATABLE ::   ztmp_fly_dta                 ! local array of values on input grid     
1350      !!----------------------------------------------------------------------
1351      ipk = SIZE(dta, 3)
1352      !
1353      !! for weighted interpolation we have weights at four corners of a box surrounding
1354      !! a model grid point, each weight is multiplied by a grid value (bilinear case)
1355      !! or by a grid value and gradients at the corner point (bicubic case)
1356      !! so we need to have a 4 by 4 subgrid surrounding each model point to cover both cases
1357
1358      !! sub grid from non-model input grid which encloses all grid points in this nemo process
1359      jpimin = ref_wgts(kw)%botleft(1)
1360      jpjmin = ref_wgts(kw)%botleft(2)
1361      jpiwid = ref_wgts(kw)%jpiwgt
1362      jpjwid = ref_wgts(kw)%jpjwgt
1363
1364      !! when reading in, expand this sub-grid by one halo point all the way round for calculating gradients
1365      rec1(1) = MAX( jpimin-1, 1 )
1366      rec1(2) = MAX( jpjmin-1, 1 )
1367      rec1(3) = 1
1368      recn(1) = MIN( jpiwid+2, ref_wgts(kw)%ddims(1)-rec1(1)+1 )
1369      recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
1370      recn(3) = kk
1371
1372      !! where we need to put it in the non-nemo grid fly_dta
1373      !! note that jpi1 and jpj1 only differ from 1 when jpimin and jpjmin are 1
1374      !! (ie at the extreme west or south of the whole input grid) and similarly for jpi2 and jpj2
1375      jpi1 = 2 + rec1(1) - jpimin
1376      jpj1 = 2 + rec1(2) - jpjmin
1377      jpi2 = jpi1 + recn(1) - 1
1378      jpj2 = jpj1 + recn(2) - 1
1379
1380
1381      IF( LEN_TRIM(lsmfile) > 0 ) THEN
1382      !! indeces for ztmp_fly_dta
1383      ! --------------------------
1384         rec1_lsm(1)=MAX(rec1(1)-nn_lsm,1)  ! starting index for enlarged external data, x direction
1385         rec1_lsm(2)=MAX(rec1(2)-nn_lsm,1)  ! starting index for enlarged external data, y direction
1386         rec1_lsm(3) = 1                    ! vertical dimension
1387         recn_lsm(1)=MIN(rec1(1)-rec1_lsm(1)+recn(1)+nn_lsm,ref_wgts(kw)%ddims(1)-rec1_lsm(1)) ! n points in x direction
1388         recn_lsm(2)=MIN(rec1(2)-rec1_lsm(2)+recn(2)+nn_lsm,ref_wgts(kw)%ddims(2)-rec1_lsm(2)) ! n points in y direction
1389         recn_lsm(3) = kk                   ! number of vertical levels in the input file
1390
1391      !  Avoid out of bound
1392         jpimin_lsm = MAX( rec1_lsm(1)+1, 1 )
1393         jpjmin_lsm = MAX( rec1_lsm(2)+1, 1 )
1394         jpiwid_lsm = MIN( recn_lsm(1)-2,ref_wgts(kw)%ddims(1)-rec1(1)+1)
1395         jpjwid_lsm = MIN( recn_lsm(2)-2,ref_wgts(kw)%ddims(2)-rec1(2)+1)
1396
1397         jpi1_lsm = 2+rec1_lsm(1)-jpimin_lsm
1398         jpj1_lsm = 2+rec1_lsm(2)-jpjmin_lsm
1399         jpi2_lsm = jpi1_lsm + recn_lsm(1) - 1
1400         jpj2_lsm = jpj1_lsm + recn_lsm(2) - 1
1401
1402
1403         itmpi=jpi2_lsm-jpi1_lsm+1
1404         itmpj=jpj2_lsm-jpj1_lsm+1
1405         itmpz=kk
1406         ALLOCATE(ztmp_fly_dta(itmpi,itmpj,itmpz))
1407         ztmp_fly_dta(:,:,:) = 0.0
1408         SELECT CASE( SIZE(ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),3) )
1409         CASE(1)
1410              CALL iom_get( num, jpdom_unknown, clvar, ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,1),   &
1411                 &          nrec, kstart = rec1_lsm, kcount = recn_lsm)
1412         CASE DEFAULT
1413              CALL iom_get( num, jpdom_unknown, clvar, ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),   &
1414                 &          nrec, kstart = rec1_lsm, kcount = recn_lsm)
1415         END SELECT
1416         CALL apply_seaoverland(lsmfile,ztmp_fly_dta(jpi1_lsm:jpi2_lsm,jpj1_lsm:jpj2_lsm,:),                  &
1417                 &                                      jpi1_lsm,jpi2_lsm,jpj1_lsm,jpj2_lsm,                  &
1418                 &                                      itmpi,itmpj,itmpz,rec1_lsm,recn_lsm)
1419
1420
1421         ! Relative indeces for remapping
1422         ii_lsm1 = (rec1(1)-rec1_lsm(1))+1
1423         ii_lsm2 = (ii_lsm1+recn(1))-1
1424         ij_lsm1 = (rec1(2)-rec1_lsm(2))+1
1425         ij_lsm2 = (ij_lsm1+recn(2))-1
1426
1427         ref_wgts(kw)%fly_dta(:,:,:) = 0.0
1428         ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:) = ztmp_fly_dta(ii_lsm1:ii_lsm2,ij_lsm1:ij_lsm2,:)
1429         DEALLOCATE(ztmp_fly_dta)
1430
1431      ELSE
1432         
1433         ref_wgts(kw)%fly_dta(:,:,:) = 0.0
1434         CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2,:), nrec, kstart = rec1, kcount = recn)
1435      ENDIF
1436     
1437
1438      !! first four weights common to both bilinear and bicubic
1439      !! data_jpi, data_jpj have already been shifted to (1,1) corresponding to botleft
1440      !! note that we have to offset by 1 into fly_dta array because of halo added to fly_dta (rec1 definition)
1441      dta(:,:,:) = 0._wp
1442      DO jn = 1,4
1443         DO_3D_00_00( 1,ipk )
1444            ni = ref_wgts(kw)%data_jpi(ji,jj,jn) + 1
1445            nj = ref_wgts(kw)%data_jpj(ji,jj,jn) + 1
1446            dta(ji,jj,jk) = dta(ji,jj,jk) + ref_wgts(kw)%data_wgt(ji,jj,jn) * ref_wgts(kw)%fly_dta(ni,nj,jk)
1447         END_3D
1448      END DO
1449
1450      IF(ref_wgts(kw)%numwgt .EQ. 16) THEN
1451
1452         !! fix up halo points that we couldnt read from file
1453         IF( jpi1 == 2 ) THEN
1454            ref_wgts(kw)%fly_dta(jpi1-1,:,:) = ref_wgts(kw)%fly_dta(jpi1,:,:)
1455         ENDIF
1456         IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
1457            ref_wgts(kw)%fly_dta(jpi2+1,:,:) = ref_wgts(kw)%fly_dta(jpi2,:,:)
1458         ENDIF
1459         IF( jpj1 == 2 ) THEN
1460            ref_wgts(kw)%fly_dta(:,jpj1-1,:) = ref_wgts(kw)%fly_dta(:,jpj1,:)
1461         ENDIF
1462         IF( jpj2 + jpjmin - 1 == ref_wgts(kw)%ddims(2)+1 .AND. jpj2 .LT. jpjwid+2 ) THEN
1463            ref_wgts(kw)%fly_dta(:,jpj2+1,:) = 2.0*ref_wgts(kw)%fly_dta(:,jpj2,:) - ref_wgts(kw)%fly_dta(:,jpj2-1,:)
1464         ENDIF
1465         
1466         !! if data grid is cyclic we can do better on east-west edges
1467         !! but have to allow for whether first and last columns are coincident
1468         IF( ref_wgts(kw)%cyclic ) THEN
1469            rec1(2) = MAX( jpjmin-1, 1 )
1470            recn(1) = 1
1471            recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
1472            jpj1 = 2 + rec1(2) - jpjmin
1473            jpj2 = jpj1 + recn(2) - 1
1474            IF( jpi1 == 2 ) THEN
1475               rec1(1) = ref_wgts(kw)%ddims(1) - ref_wgts(kw)%overlap
1476               CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, kstart = rec1, kcount = recn)
1477               ref_wgts(kw)%fly_dta(jpi1-1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
1478            ENDIF
1479            IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
1480               rec1(1) = 1 + ref_wgts(kw)%overlap
1481               CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col(:,jpj1:jpj2,:), nrec, kstart = rec1, kcount = recn)
1482               ref_wgts(kw)%fly_dta(jpi2+1,jpj1:jpj2,:) = ref_wgts(kw)%col(1,jpj1:jpj2,:)
1483            ENDIF
1484         ENDIF
1485         !
1486!!$         DO jn = 1,4
1487!!$            DO_3D_00_00( 1,ipk )
1488!!$               ni = ref_wgts(kw)%data_jpi(ji,jj,jn) + 1
1489!!$               nj = ref_wgts(kw)%data_jpj(ji,jj,jn) + 1
1490!!$               dta(ji,jj,jk) = dta(ji,jj,jk)   &
1491!!$                  ! gradient in the i direction
1492!!$                  &            + ref_wgts(kw)%data_wgt(ji,jj,jn+4) * 0.5_wp *                                    &
1493!!$                  &                (ref_wgts(kw)%fly_dta(ni+1,nj  ,jk) - ref_wgts(kw)%fly_dta(ni-1,nj  ,jk))     &
1494!!$                  ! gradient in the j direction
1495!!$                  &            + ref_wgts(kw)%data_wgt(ji,jj,jn+8) * 0.5_wp *                                    &
1496!!$                  &                (ref_wgts(kw)%fly_dta(ni  ,nj+1,jk) - ref_wgts(kw)%fly_dta(ni  ,nj-1,jk))     &
1497!!$                  ! gradient in the ij direction
1498!!$                  &            + ref_wgts(kw)%data_wgt(ji,jj,jn+12) * 0.25_wp *                                  &
1499!!$                  &               ((ref_wgts(kw)%fly_dta(ni+1,nj+1,jk) - ref_wgts(kw)%fly_dta(ni-1,nj+1,jk)) -   &
1500!!$                  &                (ref_wgts(kw)%fly_dta(ni+1,nj-1,jk) - ref_wgts(kw)%fly_dta(ni-1,nj-1,jk)))
1501!!$            END_3D
1502!!$         END DO
1503         !
1504         DO jn = 1,4
1505            DO_3D_00_00( 1,ipk )
1506               ni = ref_wgts(kw)%data_jpi(ji,jj,jn)
1507               nj = ref_wgts(kw)%data_jpj(ji,jj,jn)
1508               ! gradient in the i direction
1509               dta(ji,jj,jk) = dta(ji,jj,jk) + ref_wgts(kw)%data_wgt(ji,jj,jn+4) * 0.5_wp *         &
1510                  &                (ref_wgts(kw)%fly_dta(ni+2,nj+1,jk) - ref_wgts(kw)%fly_dta(ni  ,nj+1,jk))
1511            END_3D
1512         END DO
1513         DO jn = 1,4
1514            DO_3D_00_00( 1,ipk )
1515               ni = ref_wgts(kw)%data_jpi(ji,jj,jn)
1516               nj = ref_wgts(kw)%data_jpj(ji,jj,jn)
1517               ! gradient in the j direction
1518               dta(ji,jj,jk) = dta(ji,jj,jk) + ref_wgts(kw)%data_wgt(ji,jj,jn+8) * 0.5_wp *         &
1519                  &                (ref_wgts(kw)%fly_dta(ni+1,nj+2,jk) - ref_wgts(kw)%fly_dta(ni+1,nj  ,jk))
1520            END_3D
1521         END DO
1522         DO jn = 1,4
1523            DO_3D_00_00( 1,ipk )
1524               ni = ref_wgts(kw)%data_jpi(ji,jj,jn)
1525               nj = ref_wgts(kw)%data_jpj(ji,jj,jn)
1526               ! gradient in the ij direction
1527               dta(ji,jj,jk) = dta(ji,jj,jk) + ref_wgts(kw)%data_wgt(ji,jj,jn+12) * 0.25_wp * (     &
1528                  &                (ref_wgts(kw)%fly_dta(ni+2,nj+2,jk) - ref_wgts(kw)%fly_dta(ni  ,nj+2,jk)) -   &
1529                  &                (ref_wgts(kw)%fly_dta(ni+2,nj  ,jk) - ref_wgts(kw)%fly_dta(ni  ,nj  ,jk)))
1530            END_3D
1531         END DO
1532         !
1533      ENDIF
1534      !
1535   END SUBROUTINE fld_interp
1536
1537
1538   FUNCTION fld_filename( sdjf, kday, kmonth, kyear )
1539      !!---------------------------------------------------------------------
1540      !!                    ***  FUNCTION fld_filename ***
1541      !!
1542      !! ** Purpose :   define the filename according to a given date
1543      !!---------------------------------------------------------------------
1544      TYPE(FLD), INTENT(in) ::   sdjf         ! input field related variables
1545      INTEGER  , INTENT(in) ::   kday, kmonth, kyear
1546      !
1547      CHARACTER(len = 256) ::   clname, fld_filename
1548      !!---------------------------------------------------------------------
1549
1550     
1551      ! build the new filename if not climatological data
1552      clname=TRIM(sdjf%clrootname)
1553      !
1554      ! note that sdjf%ln_clim is is only acting on the presence of the year in the file name
1555      IF( .NOT. sdjf%ln_clim ) THEN   
1556                                         WRITE(clname, '(a,"_y",i4.4)' ) TRIM( sdjf%clrootname ), kyear    ! add year
1557         IF( sdjf%clftyp /= 'yearly' )   WRITE(clname, '(a, "m",i2.2)' ) TRIM( clname          ), kmonth   ! add month
1558      ELSE
1559         ! build the new filename if climatological data
1560         IF( sdjf%clftyp /= 'yearly' )   WRITE(clname, '(a,"_m",i2.2)' ) TRIM( sdjf%clrootname ), kmonth   ! add month
1561      ENDIF
1562      IF(    sdjf%clftyp == 'daily' .OR. sdjf%clftyp(1:4) == 'week' ) &
1563         &                               WRITE(clname, '(a,"d" ,i2.2)' ) TRIM( clname          ), kday     ! add day
1564
1565      fld_filename = clname
1566     
1567   END FUNCTION fld_filename
1568
1569
1570   FUNCTION ksec_week( cdday )
1571      !!---------------------------------------------------------------------
1572      !!                    ***  FUNCTION ksec_week ***
1573      !!
1574      !! ** Purpose :   seconds between 00h of the beginning of the week and half of the current time step
1575      !!---------------------------------------------------------------------
1576      CHARACTER(len=*), INTENT(in)   ::   cdday   ! first 3 letters of the first day of the weekly file
1577      !!
1578      INTEGER                        ::   ksec_week      ! output variable
1579      INTEGER                        ::   ijul, ishift   ! local integer
1580      CHARACTER(len=3),DIMENSION(7)  ::   cl_week 
1581      !!----------------------------------------------------------------------
1582      cl_week = (/"sun","sat","fri","thu","wed","tue","mon"/)
1583      DO ijul = 1, 7
1584         IF( cl_week(ijul) == TRIM(cdday) ) EXIT
1585      END DO
1586      IF( ijul .GT. 7 )   CALL ctl_stop( 'ksec_week: wrong day for sdjf%clftyp(6:8): '//TRIM(cdday) )
1587      !
1588      ishift = ijul * NINT(rday)
1589      !
1590      ksec_week = nsec_monday + ishift
1591      ksec_week = MOD( ksec_week, 7*NINT(rday) )
1592      !
1593   END FUNCTION ksec_week
1594
1595   !!======================================================================
1596END MODULE fldread
Note: See TracBrowser for help on using the repository browser.