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

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

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

Last change on this file since 13463 was 13463, checked in by andmirek, 4 years ago

Ticket #2195:update to trunk 13461

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