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/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/SBC – NEMO

source: NEMO/branches/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/SBC/fldread.F90 @ 12866

Last change on this file since 12866 was 12866, checked in by smasson, 4 years ago

Extra_Halo: using input files without halos, see #2366

  • Property svn:keywords set to Id
File size: 90.0 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)   ::   cltype      ! 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)              ::   cltype       ! type of data file 'daily', 'monthly' or yearly'
72      INTEGER                         ::   num          ! iom id of the jpfld files to be read
73      INTEGER , DIMENSION(2,2)        ::   nrec         ! before/after record (1: index, 2: second since Jan. 1st 00h of yr nit000)
74      INTEGER                         ::   nbb          ! index of before values
75      INTEGER                         ::   naa          ! index of after  values
76      INTEGER , ALLOCATABLE, DIMENSION(:      ) ::   nrecsec   !
77      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:  ) ::   fnow   ! input fields interpolated to now time step
78      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   fdta   ! 2 consecutive record of input fields
79      CHARACTER(len = 256)            ::   wgtname      ! current name of the NetCDF weight file acting as a key
80      !                                                 ! into the WGTLIST structure
81      CHARACTER(len = 34)             ::   vcomp        ! symbolic name for a vector component that needs rotation
82      LOGICAL, DIMENSION(2)           ::   rotn         ! flag to indicate whether before/after field has been rotated
83      INTEGER                         ::   nreclast     ! last record to be read in the current file
84      CHARACTER(len = 256)            ::   lsmname      ! current name of the NetCDF mask file acting as a key
85      !                                                 !
86      !                                                 ! Variables related to BDY
87      INTEGER                         ::   igrd         !   grid type for bdy data
88      INTEGER                         ::   ibdy         !   bdy set id number
89      INTEGER, POINTER, DIMENSION(:)  ::   imap         !   Array of integer pointers to 1D arrays
90      LOGICAL                         ::   ltotvel      !   total velocity or not (T/F)
91      LOGICAL                         ::   lzint        !   T if it requires a vertical interpolation
92   END TYPE FLD
93
94!$AGRIF_DO_NOT_TREAT
95
96   !! keep list of all weights variables so they're only read in once
97   !! need to add AGRIF directives not to process this structure
98   !! also need to force wgtname to include AGRIF nest number
99   TYPE         ::   WGT        !: Input weights related variables
100      CHARACTER(len = 256)                    ::   wgtname      ! current name of the NetCDF weight file
101      INTEGER , DIMENSION(2)                  ::   ddims        ! shape of input grid
102      INTEGER , DIMENSION(2)                  ::   botleft      ! top left corner of box in input grid containing
103      !                                                         ! current processor grid
104      INTEGER , DIMENSION(2)                  ::   topright     ! top right corner of box
105      INTEGER                                 ::   jpiwgt       ! width of box on input grid
106      INTEGER                                 ::   jpjwgt       ! height of box on input grid
107      INTEGER                                 ::   numwgt       ! number of weights (4=bilinear, 16=bicubic)
108      INTEGER                                 ::   nestid       ! for agrif, keep track of nest we're in
109      INTEGER                                 ::   overlap      ! =0 when cyclic grid has no overlapping EW columns
110      !                                                         ! =>1 when they have one or more overlapping columns     
111      !                                                         ! =-1 not cyclic
112      LOGICAL                                 ::   cyclic       ! east-west cyclic or not
113      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpi     ! array of source integers
114      INTEGER,  DIMENSION(:,:,:), POINTER     ::   data_jpj     ! array of source integers
115      REAL(wp), DIMENSION(:,:,:), POINTER     ::   data_wgt     ! array of weights on model grid
116      REAL(wp), DIMENSION(:,:,:), POINTER     ::   fly_dta      ! array of values on input grid
117      REAL(wp), DIMENSION(:,:,:), POINTER     ::   col          ! temporary array for reading in columns
118   END TYPE WGT
119
120   INTEGER,     PARAMETER             ::   tot_wgts = 20
121   TYPE( WGT ), DIMENSION(tot_wgts)   ::   ref_wgts     ! array of wgts
122   INTEGER                            ::   nxt_wgt = 1  ! point to next available space in ref_wgts array
123   INTEGER                            ::   nflag = 0
124   REAL(wp), PARAMETER                ::   undeff_lsm = -999.00_wp
125
126!$AGRIF_END_DO_NOT_TREAT
127
128   !! * Substitutions
129#  include "do_loop_substitute.h90"
130   !!----------------------------------------------------------------------
131   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
132   !! $Id$
133   !! Software governed by the CeCILL license (see ./LICENSE)
134   !!----------------------------------------------------------------------
135CONTAINS
136
137   SUBROUTINE fld_read( kt, kn_fsbc, sd, kit, pt_offset, Kmm )
138      !!---------------------------------------------------------------------
139      !!                    ***  ROUTINE fld_read  ***
140      !!                   
141      !! ** Purpose :   provide at each time step the surface ocean fluxes
142      !!                (momentum, heat, freshwater and runoff)
143      !!
144      !! ** Method  :   READ each input fields in NetCDF files using IOM
145      !!      and intepolate it to the model time-step.
146      !!         Several assumptions are made on the input file:
147      !!      blahblahblah....
148      !!----------------------------------------------------------------------
149      INTEGER  , INTENT(in   )               ::   kt        ! ocean time step
150      INTEGER  , INTENT(in   )               ::   kn_fsbc   ! sbc computation period (in time step)
151      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
152      INTEGER  , INTENT(in   ), OPTIONAL     ::   kit       ! subcycle timestep for timesplitting option
153      REAL(wp) , INTENT(in   ), OPTIONAL     ::   pt_offset ! provide fields at time other than "now"
154      INTEGER  , INTENT(in   ), OPTIONAL     ::   Kmm       ! ocean time level index
155      !!
156      INTEGER  ::   imf          ! size of the structure sd
157      INTEGER  ::   jf           ! dummy indices
158      INTEGER  ::   isecsbc      ! number of seconds between Jan. 1st 00h of nit000 year and the middle of sbc time step
159      INTEGER  ::   ibb, iaa     ! shorter name for sd(jf)%nbb and sd(jf)%naa
160      LOGICAL  ::   ll_firstcall ! true if this is the first call to fld_read for this set of fields
161      REAL(wp) ::   zt_offset    ! local time offset variable
162      REAL(wp) ::   ztinta       ! ratio applied to after  records when doing time interpolation
163      REAL(wp) ::   ztintb       ! ratio applied to before records when doing time interpolation
164      CHARACTER(LEN=1000) ::   clfmt  ! write format
165      !!---------------------------------------------------------------------
166      ll_firstcall = kt == nit000
167      IF( PRESENT(kit) )   ll_firstcall = ll_firstcall .and. kit == 1
168
169      IF( nn_components == jp_iam_sas ) THEN   ;   zt_offset = REAL( nn_fsbc, wp )
170      ELSE                                     ;   zt_offset = 0.
171      ENDIF
172      IF( PRESENT(pt_offset) )   zt_offset = pt_offset
173
174      ! Note that all varibles starting by nsec_* are shifted time by +1/2 time step to be centrered
175      IF( PRESENT(kit) ) THEN   ! ignore kn_fsbc in this case
176         isecsbc = nsec_year + nsec1jan000 + NINT( (     REAL(      kit,wp) + zt_offset ) * rn_Dt / REAL(nn_e,wp) )
177      ELSE                      ! middle of sbc time step
178         ! note: we use kn_fsbc-1 because nsec_year is defined at the middle of the current time step
179         isecsbc = nsec_year + nsec1jan000 + NINT( ( 0.5*REAL(kn_fsbc-1,wp) + zt_offset ) * rn_Dt )
180      ENDIF
181      imf = SIZE( sd )
182      !
183      IF( ll_firstcall ) THEN                      ! initialization
184         DO jf = 1, imf 
185            IF( TRIM(sd(jf)%clrootname) == 'NOT USED' )   CYCLE
186            CALL fld_init( isecsbc, sd(jf) )       ! read each before field (put them in after as they will be swapped)
187         END DO
188         IF( lwp ) CALL wgt_print()                ! control print
189      ENDIF
190      !                                            ! ====================================== !
191      IF( MOD( kt-1, kn_fsbc ) == 0 ) THEN         ! update field at each kn_fsbc time-step !
192         !                                         ! ====================================== !
193         !
194         DO jf = 1, imf                            ! ---   loop over field   --- !
195            !
196            IF( TRIM(sd(jf)%clrootname) == 'NOT USED' )   CYCLE
197            CALL fld_update( isecsbc, sd(jf), Kmm )
198            !
199         END DO                                    ! --- end loop over field --- !
200
201         CALL fld_rot( kt, sd )                    ! rotate vector before/now/after fields if needed
202
203         DO jf = 1, imf                            ! ---   loop over field   --- !
204            !
205            IF( TRIM(sd(jf)%clrootname) == 'NOT USED' )   CYCLE
206            !
207            ibb = sd(jf)%nbb   ;   iaa = sd(jf)%naa
208            !
209            IF( sd(jf)%ln_tint ) THEN              ! temporal interpolation
210               IF(lwp .AND. kt - nit000 <= 100 ) THEN
211                  clfmt = "('   fld_read: var ', a, ' kt = ', i8, ' (', f9.4,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
212                     &    "', records b/a: ', i6.4, '/', i6.4, ' (days ', f9.4,'/', f9.4, ')')"
213                  WRITE(numout, clfmt)  TRIM( sd(jf)%clvar ), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,   &           
214                     & 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
215                  WRITE(numout, *) '      zt_offset is : ',zt_offset
216               ENDIF
217               ! temporal interpolation weights
218               ztinta =  REAL( isecsbc - sd(jf)%nrec(2,ibb), wp ) / REAL( sd(jf)%nrec(2,iaa) - sd(jf)%nrec(2,ibb), wp )
219               ztintb =  1. - ztinta
220               sd(jf)%fnow(:,:,:) = ztintb * sd(jf)%fdta(:,:,:,ibb) + ztinta * sd(jf)%fdta(:,:,:,iaa)
221            ELSE   ! nothing to do...
222               IF(lwp .AND. kt - nit000 <= 100 ) THEN
223                  clfmt = "('   fld_read: var ', a, ' kt = ', i8,' (', f9.4,' days), Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
224                     &    "', record: ', i6.4, ' (days ', f9.4, ' <-> ', f9.4, ')')"
225                  WRITE(numout, clfmt) TRIM(sd(jf)%clvar), kt, REAL(isecsbc,wp)/rday, nyear, nmonth, nday,    &
226                     &                 sd(jf)%nrec(1,iaa), REAL(sd(jf)%nrec(2,ibb),wp)/rday, REAL(sd(jf)%nrec(2,iaa),wp)/rday
227               ENDIF
228            ENDIF
229            !
230            IF( kt == nitend - kn_fsbc + 1 )   CALL iom_close( sd(jf)%num )   ! Close the input files
231
232         END DO                                    ! --- end loop over field --- !
233         !
234      ENDIF
235      !
236   END SUBROUTINE fld_read
237
238
239   SUBROUTINE fld_init( ksecsbc, sdjf )
240      !!---------------------------------------------------------------------
241      !!                    ***  ROUTINE fld_init  ***
242      !!
243      !! ** Purpose :  - first call(s) to fld_def to define before values
244      !!               - open file
245      !!----------------------------------------------------------------------
246      INTEGER  , INTENT(in   ) ::   ksecsbc   !
247      TYPE(FLD), INTENT(inout) ::   sdjf         ! input field related variables
248      !!---------------------------------------------------------------------
249      !
250      IF( nflag == 0 )   nflag = -( HUGE(0) - 10 )
251      !
252      CALL fld_def( sdjf )
253      IF( sdjf%ln_tint .AND. ksecsbc < sdjf%nrecsec(1) )   CALL fld_def( sdjf, ldprev = .TRUE. )
254      !
255      CALL fld_clopn( sdjf )
256      sdjf%nrec(:,sdjf%naa) = (/ 1, nflag /)  ! default definition to force flp_update to read the file.
257      !
258   END SUBROUTINE fld_init
259
260
261   SUBROUTINE fld_update( ksecsbc, sdjf, Kmm )
262      !!---------------------------------------------------------------------
263      !!                    ***  ROUTINE fld_update  ***
264      !!
265      !! ** Purpose : Compute
266      !!              if sdjf%ln_tint = .TRUE.
267      !!                  nrec(:,iaa): record number and its time (nrec(:,ibb) is obtained from nrec(:,iaa) when swapping)
268      !!              if sdjf%ln_tint = .FALSE.
269      !!                  nrec(1,iaa): record number
270      !!                  nrec(2,ibb) and nrec(2,iaa): time of the beginning and end of the record
271      !!----------------------------------------------------------------------
272      INTEGER  ,           INTENT(in   ) ::   ksecsbc   !
273      TYPE(FLD),           INTENT(inout) ::   sdjf      ! input field related variables
274      INTEGER  , OPTIONAL, INTENT(in   ) ::   Kmm    ! ocean time level index
275      !
276      INTEGER  ::   ja           ! end of this record (in seconds)
277      INTEGER  ::   ibb, iaa     ! shorter name for sdjf%nbb and sdjf%naa
278      !!----------------------------------------------------------------------
279      ibb = sdjf%nbb   ;   iaa = sdjf%naa
280      !
281      IF( ksecsbc > sdjf%nrec(2,iaa) ) THEN     ! --> we need to update after data
282       
283         ! find where is the new after record... (it is not necessary sdjf%nrec(1,iaa)+1 )
284         ja = sdjf%nrec(1,iaa)
285         DO WHILE ( ksecsbc >= sdjf%nrecsec(ja) .AND. ja < sdjf%nreclast )   ! Warning: make sure ja <= sdjf%nreclast in this test
286            ja = ja + 1
287         END DO
288         IF( ksecsbc > sdjf%nrecsec(ja) )   ja = ja + 1   ! in case ksecsbc > sdjf%nrecsec(sdjf%nreclast)
289
290         ! if ln_tint and if the new after is not ja+1, we need also to update after data before the swap
291         ! so, after the swap, sdjf%nrec(2,ibb) will still be the closest value located just before ksecsbc
292         IF( sdjf%ln_tint .AND. ( ja > sdjf%nrec(1,iaa) + 1 .OR. sdjf%nrec(2,iaa) == nflag ) ) THEN
293            sdjf%nrec(:,iaa) = (/ ja-1, sdjf%nrecsec(ja-1) /)   ! update nrec(:,iaa) with before information
294            CALL fld_get( sdjf, Kmm )                           ! read after data that will be used as before data
295         ENDIF
296           
297         ! if after is in the next file...
298         IF( ja > sdjf%nreclast ) THEN
299           
300            CALL fld_def( sdjf )
301            IF( ksecsbc > sdjf%nrecsec(sdjf%nreclast) )   CALL fld_def( sdjf, ldnext = .TRUE. )
302            CALL fld_clopn( sdjf )           ! open next file
303           
304            ! find where is after in this new file
305            ja = 1
306            DO WHILE ( ksecsbc > sdjf%nrecsec(ja) .AND. ja < sdjf%nreclast )
307               ja = ja + 1
308            END DO
309            IF( ksecsbc > sdjf%nrecsec(ja) )   ja = ja + 1   ! in case ksecsbc > sdjf%nrecsec(sdjf%nreclast)
310           
311            IF( ja > sdjf%nreclast ) THEN
312               CALL ctl_stop( "STOP", "fld_def: need next-next file? we should not be there... file: "//TRIM(sdjf%clrootname) )
313            ENDIF
314           
315            ! if ln_tint and if after is not the first record, we must (potentially again) update after data before the swap
316            IF( sdjf%ln_tint .AND. ja > 1 ) THEN
317               IF( sdjf%nrecsec(0) /= nflag ) THEN                    ! no trick used: after file is not the current file
318                  sdjf%nrec(:,iaa) = (/ ja-1, sdjf%nrecsec(ja-1) /)   ! update nrec(:,iaa) with before information
319                  CALL fld_get( sdjf, Kmm )                           ! read after data that will be used as before data
320               ENDIF
321            ENDIF
322           
323         ENDIF
324
325         IF( sdjf%ln_tint ) THEN                                ! Swap data
326            sdjf%nbb = sdjf%naa                                 !    swap indices
327            sdjf%naa = 3 - sdjf%naa                             !    = 2(1) if naa == 1(2)
328         ELSE                                                   ! No swap
329            sdjf%nrec(:,ibb) = (/ ja-1, sdjf%nrecsec(ja-1) /)   !    only for print
330         ENDIF
331           
332         ! read new after data
333         sdjf%nrec(:,sdjf%naa) = (/ ja, sdjf%nrecsec(ja) /)     ! update nrec(:,naa) as it is used by fld_get
334         CALL fld_get( sdjf, Kmm )                              ! read after data (with nrec(:,naa) informations)
335       
336      ENDIF
337      !
338   END SUBROUTINE fld_update
339
340
341   SUBROUTINE fld_get( sdjf, Kmm )
342      !!---------------------------------------------------------------------
343      !!                    ***  ROUTINE fld_get  ***
344      !!
345      !! ** Purpose :   read the data
346      !!----------------------------------------------------------------------
347      TYPE(FLD),           INTENT(inout) ::   sdjf   ! input field related variables
348      INTEGER  , OPTIONAL, INTENT(in   ) ::   Kmm    ! ocean time level index
349      !
350      INTEGER ::   ipk      ! number of vertical levels of sdjf%fdta ( 2D: ipk=1 ; 3D: ipk=jpk )
351      INTEGER ::   iaa      ! shorter name for sdjf%naa
352      INTEGER ::   iw       ! index into wgts array
353      INTEGER ::   ipdom    ! index of the domain
354      INTEGER ::   idvar    ! variable ID
355      INTEGER ::   idmspc   ! number of spatial dimensions
356      LOGICAL ::   lmoor    ! C1D case: point data
357      !!---------------------------------------------------------------------
358      !
359      ipk = SIZE( sdjf%fnow, 3 )
360      iaa = sdjf%naa
361      !
362      IF( ASSOCIATED(sdjf%imap) ) THEN
363         IF( sdjf%ln_tint ) THEN   ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fdta(:,:,:,iaa), sdjf%nrec(1,iaa),   &
364            &                                        sdjf%imap, sdjf%igrd, sdjf%ibdy, sdjf%ltotvel, sdjf%lzint, Kmm )
365         ELSE                      ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fnow(:,:,:    ), sdjf%nrec(1,iaa),   &
366            &                                        sdjf%imap, sdjf%igrd, sdjf%ibdy, sdjf%ltotvel, sdjf%lzint, Kmm )
367         ENDIF
368      ELSE IF( LEN(TRIM(sdjf%wgtname)) > 0 ) THEN
369         CALL wgt_list( sdjf, iw )
370         IF( sdjf%ln_tint ) THEN   ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw, ipk, sdjf%fdta(:,:,:,iaa),   & 
371            &                                                                          sdjf%nrec(1,iaa), sdjf%lsmname )
372            CALL lbc_lnk( 'fldread', sdjf%fdta(:,:,:,iaa), 'T', 1._wp, kfillmode = jpfillcopy )
373         ELSE                      ;   CALL fld_interp( sdjf%num, sdjf%clvar, iw, ipk, sdjf%fnow(:,:,:    ),   &
374            &                                                                          sdjf%nrec(1,iaa), sdjf%lsmname )
375            CALL lbc_lnk( 'fldread', sdjf%fnow(:,:,:    ), 'T', 1._wp, kfillmode = jpfillcopy )
376         ENDIF
377      ELSE
378         IF( SIZE(sdjf%fnow, 1) == jpi ) THEN   ;   ipdom = jpdom_global
379         ELSE                                   ;   ipdom = jpdom_unknown
380         ENDIF
381         ! C1D case: If product of spatial dimensions == ipk, then x,y are of
382         ! size 1 (point/mooring data): this must be read onto the central grid point
383         idvar  = iom_varid( sdjf%num, sdjf%clvar )
384         idmspc = iom_file ( sdjf%num )%ndims( idvar )
385         IF( iom_file( sdjf%num )%luld( idvar ) )   idmspc = idmspc - 1
386         lmoor  = (  idmspc == 0 .OR. PRODUCT( iom_file( sdjf%num )%dimsz( 1:MAX(idmspc,1) ,idvar ) ) == ipk  )
387         !
388         SELECT CASE( ipk )
389         CASE(1)
390            IF( lk_c1d .AND. lmoor ) THEN
391               IF( sdjf%ln_tint ) THEN
392                  CALL iom_get( sdjf%num, sdjf%clvar, sdjf%fdta(2,2,1,iaa), sdjf%nrec(1,iaa) )
393                  CALL lbc_lnk( 'fldread', sdjf%fdta(:,:,1,iaa),'T',1., kfillmode = jpfillcopy )
394               ELSE
395                  CALL iom_get( sdjf%num, sdjf%clvar, sdjf%fnow(2,2,1    ), sdjf%nrec(1,iaa) )
396                  CALL lbc_lnk( 'fldread', sdjf%fnow(:,:,1    ),'T',1., kfillmode = jpfillcopy )
397               ENDIF
398            ELSE
399               IF( sdjf%ln_tint ) THEN
400                  CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fdta(:,:,1,iaa), sdjf%nrec(1,iaa), kfill = jpfillcopy )
401               ELSE
402                  CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fnow(:,:,1    ), sdjf%nrec(1,iaa), kfill = jpfillcopy )
403               ENDIF
404            ENDIF
405         CASE DEFAULT
406            IF(lk_c1d .AND. lmoor ) THEN
407               IF( sdjf%ln_tint ) THEN
408                  CALL iom_get( sdjf%num, jpdom_unknown, sdjf%clvar, sdjf%fdta(2,2,:,iaa), sdjf%nrec(1,iaa) )
409                  CALL lbc_lnk( 'fldread', sdjf%fdta(:,:,:,iaa),'T',1., kfillmode = jpfillcopy )
410               ELSE
411                  CALL iom_get( sdjf%num, jpdom_unknown, sdjf%clvar, sdjf%fnow(2,2,:   ), sdjf%nrec(1,iaa) )
412                  CALL lbc_lnk( 'fldread', sdjf%fnow(:,:,:   ),'T',1., kfillmode = jpfillcopy )
413               ENDIF
414            ELSE
415               IF( sdjf%ln_tint ) THEN
416                  CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fdta(:,:,:,iaa), sdjf%nrec(1,iaa), kfill = jpfillcopy )
417               ELSE
418                  CALL iom_get( sdjf%num, ipdom, sdjf%clvar, sdjf%fnow(:,:,:    ), sdjf%nrec(1,iaa), kfill = jpfillcopy )
419               ENDIF
420            ENDIF
421         END SELECT
422      ENDIF
423      !
424      sdjf%rotn(iaa) = .false.   ! vector not yet rotated
425      !
426   END SUBROUTINE fld_get
427
428   
429   SUBROUTINE fld_map( knum, cdvar, pdta, krec, kmap, kgrd, kbdy, ldtotvel, ldzint, Kmm )
430      !!---------------------------------------------------------------------
431      !!                    ***  ROUTINE fld_map  ***
432      !!
433      !! ** Purpose :   read global data from file and map onto local data
434      !!                using a general mapping (for open boundaries)
435      !!----------------------------------------------------------------------
436      INTEGER                   , INTENT(in   ) ::   knum         ! stream number
437      CHARACTER(LEN=*)          , INTENT(in   ) ::   cdvar        ! variable name
438      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pdta         ! bdy output field on model grid
439      INTEGER                   , INTENT(in   ) ::   krec         ! record number to read (ie time slice)
440      INTEGER , DIMENSION(:)    , INTENT(in   ) ::   kmap         ! global-to-local bdy mapping indices
441      ! optional variables used for vertical interpolation:
442      INTEGER, OPTIONAL         , INTENT(in   ) ::   kgrd         ! grid type (t, u, v)
443      INTEGER, OPTIONAL         , INTENT(in   ) ::   kbdy         ! bdy number
444      LOGICAL, OPTIONAL         , INTENT(in   ) ::   ldtotvel     ! true if total ( = barotrop + barocline) velocity
445      LOGICAL, OPTIONAL         , INTENT(in   ) ::   ldzint       ! true if 3D variable requires a vertical interpolation
446      INTEGER, OPTIONAL         , INTENT(in   ) ::   Kmm          ! ocean time level index
447      !!
448      INTEGER                                   ::   ipi          ! length of boundary data on local process
449      INTEGER                                   ::   ipj          ! length of dummy dimension ( = 1 )
450      INTEGER                                   ::   ipk          ! number of vertical levels of pdta ( 2D: ipk=1 ; 3D: ipk=jpk )
451      INTEGER                                   ::   ipkb         ! number of vertical levels in boundary data file
452      INTEGER                                   ::   idvar        ! variable ID
453      INTEGER                                   ::   indims       ! number of dimensions of the variable
454      INTEGER, DIMENSION(4)                     ::   idimsz       ! size of variable dimensions
455      REAL(wp)                                  ::   zfv          ! fillvalue
456      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zz_read      ! work space for global boundary data
457      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zdta_read    ! work space local data requiring vertical interpolation
458      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zdta_read_z  ! work space local data requiring vertical interpolation
459      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::   zdta_read_dz ! work space local data requiring vertical interpolation
460      CHARACTER(LEN=1),DIMENSION(3)             ::   clgrid
461      LOGICAL                                   ::   lluld        ! is the variable using the unlimited dimension
462      LOGICAL                                   ::   llzint       ! local value of ldzint
463      !!---------------------------------------------------------------------
464      !
465      clgrid = (/'t','u','v'/)
466      !
467      ipi = SIZE( pdta, 1 )
468      ipj = SIZE( pdta, 2 )   ! must be equal to 1
469      ipk = SIZE( pdta, 3 )
470      !
471      llzint = .FALSE.
472      IF( PRESENT(ldzint) )   llzint = ldzint
473      !
474      idvar = iom_varid( knum, cdvar, kndims = indims, kdimsz = idimsz, lduld = lluld  )
475      IF( indims == 4 .OR. ( indims == 3 .AND. .NOT. lluld ) ) THEN   ;   ipkb = idimsz(3)   ! xy(zl)t or xy(zl)
476      ELSE                                                            ;   ipkb = 1           ! xy or xyt
477      ENDIF
478      !
479      ALLOCATE( zz_read( idimsz(1), idimsz(2), ipkb ) )  ! ++++++++ !!! this can be very big...         
480      !
481      IF( ipk == 1 ) THEN
482
483         IF( ipkb /= 1 ) CALL ctl_stop( 'fld_map : we must have ipkb = 1 to read surface data' )
484         CALL iom_get ( knum, jpdom_unknown, cdvar, zz_read(:,:,1), krec )   ! call iom_get with a 2D file
485         CALL fld_map_core( zz_read, kmap, pdta )
486
487      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
488      ! Do we include something here to adjust barotropic velocities !
489      ! in case of a depth difference between bdy files and          !
490      ! bathymetry in the case ln_totvel = .false. and ipkb>0?       !
491      ! [as the enveloping and parital cells could change H]         !
492      !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
493
494      ELSE
495         !
496         CALL iom_get ( knum, jpdom_unknown, cdvar, zz_read(:,:,:), krec )   ! call iom_get with a 3D file
497         !
498         IF( ipkb /= ipk .OR. llzint ) THEN   ! boundary data not on model vertical grid : vertical interpolation
499            !
500            IF( ipk == jpk .AND. iom_varid(knum,'gdep'//clgrid(kgrd)) /= -1 .AND. iom_varid(knum,'e3'//clgrid(kgrd)) /= -1 ) THEN
501               
502               ALLOCATE( zdta_read(ipi,ipj,ipkb), zdta_read_z(ipi,ipj,ipkb), zdta_read_dz(ipi,ipj,ipkb) )
503               
504               CALL fld_map_core( zz_read, kmap, zdta_read )
505               CALL iom_get ( knum, jpdom_unknown, 'gdep'//clgrid(kgrd), zz_read )   ! read only once? Potential temporal evolution?
506               CALL fld_map_core( zz_read, kmap, zdta_read_z )
507               CALL iom_get ( knum, jpdom_unknown,   'e3'//clgrid(kgrd), zz_read )   ! read only once? Potential temporal evolution?
508               CALL fld_map_core( zz_read, kmap, zdta_read_dz )
509               
510               CALL iom_getatt(knum, '_FillValue', zfv, cdvar=cdvar )
511               CALL fld_bdy_interp(zdta_read, zdta_read_z, zdta_read_dz, pdta, kgrd, kbdy, zfv, ldtotvel, Kmm)
512               DEALLOCATE( zdta_read, zdta_read_z, zdta_read_dz )
513               
514            ELSE
515               IF( ipk /= jpk ) CALL ctl_stop( 'fld_map : this should be an impossible case...' )
516               WRITE(ctmp1,*) 'fld_map : vertical interpolation for bdy variable '//TRIM(cdvar)//' requires ' 
517               IF( iom_varid(knum, 'gdep'//clgrid(kgrd)) == -1 ) CALL ctl_stop( ctmp1//'gdep'//clgrid(kgrd)//' variable' )
518               IF( iom_varid(knum,   'e3'//clgrid(kgrd)) == -1 ) CALL ctl_stop( ctmp1//  'e3'//clgrid(kgrd)//' variable' )
519
520            ENDIF
521            !
522         ELSE                            ! bdy data assumed to be the same levels as bdy variables
523            !
524            CALL fld_map_core( zz_read, kmap, pdta )
525            !
526         ENDIF   ! ipkb /= ipk
527      ENDIF   ! ipk == 1
528     
529      DEALLOCATE( zz_read )
530
531   END SUBROUTINE fld_map
532
533     
534   SUBROUTINE fld_map_core( pdta_read, kmap, pdta_bdy )
535      !!---------------------------------------------------------------------
536      !!                    ***  ROUTINE fld_map_core  ***
537      !!
538      !! ** Purpose :  inner core of fld_map
539      !!----------------------------------------------------------------------
540      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read    ! global boundary data
541      INTEGER,  DIMENSION(:    ), INTENT(in   ) ::   kmap         ! global-to-local bdy mapping indices
542      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pdta_bdy     ! bdy output field on model grid
543      !!
544      INTEGER,  DIMENSION(3) ::   idim_read,  idim_bdy            ! arrays dimensions
545      INTEGER                ::   ji, jj, jk, jb                  ! loop counters
546      INTEGER                ::   im1
547      !!---------------------------------------------------------------------
548      !
549      idim_read = SHAPE( pdta_read )
550      idim_bdy  = SHAPE( pdta_bdy  )
551      !
552      ! in all cases: idim_bdy(2) == 1 .AND. idim_read(1) * idim_read(2) == idim_bdy(1)
553      ! structured BDY with rimwidth > 1                     : idim_read(2) == rimwidth /= 1
554      ! structured BDY with rimwidth == 1 or unstructured BDY: idim_read(2) == 1
555      !
556      IF( idim_read(2) > 1 ) THEN    ! structured BDY with rimwidth > 1 
557         DO jk = 1, idim_bdy(3)
558            DO jb = 1, idim_bdy(1)
559               im1 = kmap(jb) - 1
560               jj = im1 / idim_read(1) + 1
561               ji = MOD( im1, idim_read(1) ) + 1
562               pdta_bdy(jb,1,jk) =  pdta_read(ji,jj,jk)
563            END DO
564         END DO
565      ELSE
566         DO jk = 1, idim_bdy(3)
567            DO jb = 1, idim_bdy(1)   ! horizontal remap of bdy data on the local bdy
568               pdta_bdy(jb,1,jk) = pdta_read(kmap(jb),1,jk)
569            END DO
570         END DO
571      ENDIF
572     
573   END SUBROUTINE fld_map_core
574   
575   SUBROUTINE fld_bdy_interp(pdta_read, pdta_read_z, pdta_read_dz, pdta, kgrd, kbdy, pfv, ldtotvel, Kmm )
576      !!---------------------------------------------------------------------
577      !!                    ***  ROUTINE fld_bdy_interp  ***
578      !!
579      !! ** Purpose :   on the fly vertical interpolation to allow the use of
580      !!                boundary data from non-native vertical grid
581      !!----------------------------------------------------------------------
582      USE bdy_oce, ONLY:  idx_bdy         ! indexing for map <-> ij transformation
583
584      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read       ! data read in bdy file
585      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read_z     ! depth of the data read in bdy file
586      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pdta_read_dz    ! thickness of the levels in bdy file
587      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   pdta            ! output field on model grid (2 dimensional)
588      REAL(wp)                  , INTENT(in   ) ::   pfv             ! fillvalue of the data read in bdy file
589      LOGICAL                   , INTENT(in   ) ::   ldtotvel        ! true if toal ( = barotrop + barocline) velocity
590      INTEGER                   , INTENT(in   ) ::   kgrd            ! grid type (t, u, v)
591      INTEGER                   , INTENT(in   ) ::   kbdy            ! bdy number
592      INTEGER, OPTIONAL         , INTENT(in   ) ::   Kmm             ! ocean time level index
593      !!
594      INTEGER                  ::   ipi                 ! length of boundary data on local process
595      INTEGER                  ::   ipkb                ! number of vertical levels in boundary data file
596      INTEGER                  ::   ipkmax              ! number of vertical levels in boundary data file where no mask
597      INTEGER                  ::   jb, ji, jj, jk, jkb ! loop counters
598      REAL(wp)                 ::   zcoef, zi           !
599      REAL(wp)                 ::   ztrans, ztrans_new  ! transports
600      REAL(wp), DIMENSION(jpk) ::   zdepth, zdhalf      ! level and half-level depth
601      !!---------------------------------------------------------------------
602     
603      ipi  = SIZE( pdta, 1 )
604      ipkb = SIZE( pdta_read, 3 )
605     
606      DO jb = 1, ipi
607         ji = idx_bdy(kbdy)%nbi(jb,kgrd)
608         jj = idx_bdy(kbdy)%nbj(jb,kgrd)
609         !
610         ! --- max jk where input data /= FillValue --- !
611         ipkmax = 1
612         DO jkb = 2, ipkb
613            IF( pdta_read(jb,1,jkb) /= pfv )   ipkmax = MAX( ipkmax, jkb )
614         END DO
615         !
616         ! --- calculate depth at t,u,v points --- !
617         SELECT CASE( kgrd )                         
618         CASE(1)            ! depth of T points:
619            zdepth(:) = gdept(ji,jj,:,Kmm)
620         CASE(2)            ! depth of U points: we must not use gdept_n as we don't want to do a communication
621            !                 --> copy what is done for gdept_n in domvvl...
622            zdhalf(1) = 0.0_wp
623            zdepth(1) = 0.5_wp * e3uw(ji,jj,1,Kmm)
624            DO jk = 2, jpk                               ! vertical sum
625               !    zcoef = umask - wumask    ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
626               !                              ! 1 everywhere from mbkt to mikt + 1 or 1 (if no isf)
627               !                              ! 0.5 where jk = mikt     
628               !!gm ???????   BUG ?  gdept_n as well as gde3w_n  does not include the thickness of ISF ??
629               zcoef = ( umask(ji,jj,jk) - wumask(ji,jj,jk) )
630               zdhalf(jk) = zdhalf(jk-1) + e3u(ji,jj,jk-1,Kmm)
631               zdepth(jk) =          zcoef  * ( zdhalf(jk  ) + 0.5_wp * e3uw(ji,jj,jk,Kmm))  &
632                  &         + (1._wp-zcoef) * ( zdepth(jk-1) +          e3uw(ji,jj,jk,Kmm))
633            END DO
634         CASE(3)            ! depth of V points: we must not use gdept_n as we don't want to do a communication
635            !                 --> copy what is done for gdept_n in domvvl...
636            zdhalf(1) = 0.0_wp
637            zdepth(1) = 0.5_wp * e3vw(ji,jj,1,Kmm)
638            DO jk = 2, jpk                               ! vertical sum
639               !    zcoef = vmask - wvmask    ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
640               !                              ! 1 everywhere from mbkt to mikt + 1 or 1 (if no isf)
641               !                              ! 0.5 where jk = mikt     
642               !!gm ???????   BUG ?  gdept_n as well as gde3w_n  does not include the thickness of ISF ??
643               zcoef = ( vmask(ji,jj,jk) - wvmask(ji,jj,jk) )
644               zdhalf(jk) = zdhalf(jk-1) + e3v(ji,jj,jk-1,Kmm)
645               zdepth(jk) =          zcoef  * ( zdhalf(jk  ) + 0.5_wp * e3vw(ji,jj,jk,Kmm))  &
646                  &         + (1._wp-zcoef) * ( zdepth(jk-1) +          e3vw(ji,jj,jk,Kmm))
647            END DO
648         END SELECT
649         !         
650         ! --- interpolate bdy data on the model grid --- !
651         DO jk = 1, jpk
652            IF(     zdepth(jk) <= pdta_read_z(jb,1,1)      ) THEN   ! above the first level of external data
653               pdta(jb,1,jk) = pdta_read(jb,1,1)
654            ELSEIF( zdepth(jk) >  pdta_read_z(jb,1,ipkmax) ) THEN   ! below the last level of external data /= FillValue
655               pdta(jb,1,jk) = pdta_read(jb,1,ipkmax)
656            ELSE                                                    ! inbetween: vertical interpolation between jkb & jkb+1
657               DO jkb = 1, ipkmax-1
658                  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
659                     zi = ( zdepth(jk) - pdta_read_z(jb,1,jkb) ) / ( pdta_read_z(jb,1,jkb+1) - pdta_read_z(jb,1,jkb) )
660                     pdta(jb,1,jk) = pdta_read(jb,1,jkb) + zi * ( pdta_read(jb,1,jkb+1) - pdta_read(jb,1,jkb) )
661                  ENDIF
662               END DO
663            ENDIF
664         END DO   ! jpk
665         !
666      END DO   ! ipi
667
668      ! --- mask data and adjust transport --- !
669      SELECT CASE( kgrd )                         
670
671      CASE(1)                                 ! mask data (probably unecessary)
672         DO jb = 1, ipi
673            ji = idx_bdy(kbdy)%nbi(jb,kgrd)
674            jj = idx_bdy(kbdy)%nbj(jb,kgrd)
675            DO jk = 1, jpk                     
676               pdta(jb,1,jk) = pdta(jb,1,jk) * tmask(ji,jj,jk)
677            END DO
678         END DO
679         
680      CASE(2)                                 ! adjust the U-transport term
681         DO jb = 1, ipi
682            ji = idx_bdy(kbdy)%nbi(jb,kgrd)
683            jj = idx_bdy(kbdy)%nbj(jb,kgrd)
684            ztrans = 0._wp
685            DO jkb = 1, ipkb                              ! calculate transport on input grid
686               IF( pdta_read(jb,1,jkb) /= pfv )   ztrans = ztrans + pdta_read(jb,1,jkb) * pdta_read_dz(jb,1,jkb)
687            ENDDO
688            ztrans_new = 0._wp
689            DO jk = 1, jpk                                ! calculate transport on model grid
690               ztrans_new = ztrans_new +      pdta(jb,1,jk ) * e3u(ji,jj,jk,Kmm ) * umask(ji,jj,jk)
691            ENDDO
692            DO jk = 1, jpk                                ! make transport correction
693               IF(ldtotvel) THEN ! bdy data are total velocity so adjust bt transport term to match input data
694                  pdta(jb,1,jk) = ( pdta(jb,1,jk) + ( ztrans - ztrans_new ) * r1_hu(ji,jj,Kmm) ) * umask(ji,jj,jk)
695               ELSE              ! we're just dealing with bc velocity so bt transport term should sum to zero
696                  pdta(jb,1,jk) =   pdta(jb,1,jk) + (  0._wp - ztrans_new ) * r1_hu(ji,jj,Kmm)   * umask(ji,jj,jk)
697               ENDIF
698            ENDDO
699         ENDDO
700
701      CASE(3)                                 ! adjust the V-transport term
702         DO jb = 1, ipi
703            ji = idx_bdy(kbdy)%nbi(jb,kgrd)
704            jj = idx_bdy(kbdy)%nbj(jb,kgrd)
705            ztrans = 0._wp
706            DO jkb = 1, ipkb                              ! calculate transport on input grid
707               IF( pdta_read(jb,1,jkb) /= pfv )   ztrans = ztrans + pdta_read(jb,1,jkb) * pdta_read_dz(jb,1,jkb)
708            ENDDO
709            ztrans_new = 0._wp
710            DO jk = 1, jpk                                ! calculate transport on model grid
711               ztrans_new = ztrans_new +      pdta(jb,1,jk ) * e3v(ji,jj,jk,Kmm ) * vmask(ji,jj,jk)
712            ENDDO
713            DO jk = 1, jpk                                ! make transport correction
714               IF(ldtotvel) THEN ! bdy data are total velocity so adjust bt transport term to match input data
715                  pdta(jb,1,jk) = ( pdta(jb,1,jk) + ( ztrans - ztrans_new ) * r1_hv(ji,jj,Kmm) ) * vmask(ji,jj,jk)
716               ELSE              ! we're just dealing with bc velocity so bt transport term should sum to zero
717                  pdta(jb,1,jk) =   pdta(jb,1,jk) + (  0._wp - ztrans_new ) * r1_hv(ji,jj,Kmm)   * vmask(ji,jj,jk)
718               ENDIF
719            ENDDO
720         ENDDO
721      END SELECT
722     
723   END SUBROUTINE fld_bdy_interp
724
725
726   SUBROUTINE fld_rot( kt, sd )
727      !!---------------------------------------------------------------------
728      !!                    ***  ROUTINE fld_rot  ***
729      !!
730      !! ** Purpose :   Vector fields may need to be rotated onto the local grid direction
731      !!----------------------------------------------------------------------
732      INTEGER                , INTENT(in   ) ::   kt   ! ocean time step
733      TYPE(FLD), DIMENSION(:), INTENT(inout) ::   sd   ! input field related variables
734      !
735      INTEGER ::   ju, jv, jk, jn  ! loop indices
736      INTEGER ::   imf             ! size of the structure sd
737      INTEGER ::   ill             ! character length
738      INTEGER ::   iv              ! indice of V component
739      CHARACTER (LEN=100)          ::   clcomp       ! dummy weight name
740      REAL(wp), DIMENSION(jpi,jpj) ::   utmp, vtmp   ! temporary arrays for vector rotation
741      !!---------------------------------------------------------------------
742      !
743      !! (sga: following code should be modified so that pairs arent searched for each time
744      !
745      imf = SIZE( sd )
746      DO ju = 1, imf
747         IF( TRIM(sd(ju)%clrootname) == 'NOT USED' )   CYCLE
748         ill = LEN_TRIM( sd(ju)%vcomp )
749         DO jn = 2-COUNT((/sd(ju)%ln_tint/)), 2
750            IF( ill > 0 .AND. .NOT. sd(ju)%rotn(jn) ) THEN   ! find vector rotations required             
751               IF( sd(ju)%vcomp(1:1) == 'U' ) THEN      ! east-west component has symbolic name starting with 'U'
752                  ! look for the north-south component which has same symbolic name but with 'U' replaced with 'V'
753                  clcomp = 'V' // sd(ju)%vcomp(2:ill)   ! works even if ill == 1
754                  iv = -1
755                  DO jv = 1, imf
756                     IF( TRIM(sd(jv)%clrootname) == 'NOT USED' )   CYCLE
757                     IF( TRIM(sd(jv)%vcomp) == TRIM(clcomp) )   iv = jv
758                  END DO
759                  IF( iv > 0 ) THEN   ! fields ju and iv are two components which need to be rotated together
760                     DO jk = 1, SIZE( sd(ju)%fnow, 3 )
761                        IF( sd(ju)%ln_tint )THEN
762                           CALL rot_rep( sd(ju)%fdta(:,:,jk,jn), sd(iv)%fdta(:,:,jk,jn), 'T', 'en->i', utmp(:,:) )
763                           CALL rot_rep( sd(ju)%fdta(:,:,jk,jn), sd(iv)%fdta(:,:,jk,jn), 'T', 'en->j', vtmp(:,:) )
764                           sd(ju)%fdta(:,:,jk,jn) = utmp(:,:)   ;   sd(iv)%fdta(:,:,jk,jn) = vtmp(:,:)
765                        ELSE
766                           CALL rot_rep( sd(ju)%fnow(:,:,jk   ), sd(iv)%fnow(:,:,jk   ), 'T', 'en->i', utmp(:,:) )
767                           CALL rot_rep( sd(ju)%fnow(:,:,jk   ), sd(iv)%fnow(:,:,jk   ), 'T', 'en->j', vtmp(:,:) )
768                           sd(ju)%fnow(:,:,jk   ) = utmp(:,:)   ;   sd(iv)%fnow(:,:,jk   ) = vtmp(:,:)
769                        ENDIF
770                     END DO
771                     sd(ju)%rotn(jn) = .TRUE.               ! vector was rotated
772                     IF( lwp .AND. kt == nit000 )   WRITE(numout,*)   &
773                        &   'fld_read: vector pair ('//TRIM(sd(ju)%clvar)//', '//TRIM(sd(iv)%clvar)//') rotated on to model grid'
774                  ENDIF
775               ENDIF
776            ENDIF
777         END DO
778       END DO
779      !
780   END SUBROUTINE fld_rot
781
782
783   SUBROUTINE fld_def( sdjf, ldprev, ldnext )
784      !!---------------------------------------------------------------------
785      !!                    ***  ROUTINE fld_def  ***
786      !!
787      !! ** Purpose :   define the record(s) of the file and its name
788      !!----------------------------------------------------------------------
789      TYPE(FLD)        , INTENT(inout) ::   sdjf       ! input field related variables
790      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldprev     !
791      LOGICAL, OPTIONAL, INTENT(in   ) ::   ldnext     !
792      !
793      INTEGER  :: jt
794      INTEGER  :: idaysec               ! number of seconds in 1 day = NINT(rday)
795      INTEGER  :: iyr, imt, idy, isecwk
796      INTEGER  :: indexyr, indexmt
797      INTEGER  :: ireclast
798      INTEGER  :: ishift, istart
799      INTEGER, DIMENSION(2)  :: isave
800      REAL(wp) :: zfreqs
801      LOGICAL  :: llprev, llnext, llstop
802      LOGICAL  :: llprevmt, llprevyr
803      LOGICAL  :: llnextmt, llnextyr
804      !!----------------------------------------------------------------------
805      idaysec = NINT(rday)
806      !
807      IF( PRESENT(ldprev) ) THEN   ;   llprev = ldprev
808      ELSE                         ;   llprev = .FALSE.
809      ENDIF
810      IF( PRESENT(ldnext) ) THEN   ;   llnext = ldnext
811      ELSE                         ;   llnext = .FALSE.
812      ENDIF
813
814      ! current file parameters
815      IF( sdjf%cltype(1:4) == 'week' ) THEN          ! find the day of the beginning of the current week
816         isecwk = ksec_week( sdjf%cltype(6:8) )     ! seconds between the beginning of the week and half of current time step
817         llprevmt = isecwk > nsec_month               ! longer time since beginning of the current week than the current month
818         llprevyr = llprevmt .AND. nmonth == 1
819         iyr = nyear  - COUNT((/llprevyr/))
820         imt = nmonth - COUNT((/llprevmt/)) + 12 * COUNT((/llprevyr/))
821         idy = nday + nmonth_len(nmonth-1) * COUNT((/llprevmt/)) - isecwk / idaysec
822         isecwk = nsec_year - isecwk              ! seconds between 00h jan 1st of current year and current week beginning
823      ELSE
824         iyr = nyear
825         imt = nmonth
826         idy = nday
827         isecwk  = 0
828      ENDIF
829
830      ! previous file parameters
831      IF( llprev ) THEN
832         IF( sdjf%cltype(1:4) == 'week'    ) THEN     ! find the day of the beginning of previous week
833            isecwk = isecwk + 7 * idaysec         ! seconds between the beginning of previous week and half of the time step
834            llprevmt = isecwk > nsec_month            ! longer time since beginning of the previous week than the current month
835            llprevyr = llprevmt .AND. nmonth == 1
836            iyr = nyear  - COUNT((/llprevyr/))
837            imt = nmonth - COUNT((/llprevmt/)) + 12 * COUNT((/llprevyr/))
838            idy = nday + nmonth_len(nmonth-1) * COUNT((/llprevmt/)) - isecwk / idaysec
839            isecwk = nsec_year - isecwk           ! seconds between 00h jan 1st of current year and previous week beginning
840         ELSE
841            idy = nday   - COUNT((/ sdjf%cltype == 'daily'                 /))
842            imt = nmonth - COUNT((/ sdjf%cltype == 'monthly' .OR. idy == 0 /))
843            iyr = nyear  - COUNT((/ sdjf%cltype == 'yearly'  .OR. imt == 0 /))
844            IF( idy == 0 ) idy = nmonth_len(imt)
845            IF( imt == 0 ) imt = 12
846            isecwk = 0
847         ENDIF
848      ENDIF
849
850      ! next file parameters
851      IF( llnext ) THEN
852         IF( sdjf%cltype(1:4) == 'week'    ) THEN     ! find the day of the beginning of next week
853            isecwk = 7 * idaysec - isecwk         ! seconds between half of the time step and the beginning of next week
854            llnextmt = isecwk > ( nmonth_len(nmonth)*idaysec - nsec_month )   ! larger than the seconds to the end of the month
855            llnextyr = llnextmt .AND. nmonth == 12
856            iyr = nyear  + COUNT((/llnextyr/))
857            imt = nmonth + COUNT((/llnextmt/)) - 12 * COUNT((/llnextyr/))
858            idy = nday - nmonth_len(nmonth) * COUNT((/llnextmt/)) + isecwk / idaysec + 1
859            isecwk = nsec_year + isecwk           ! seconds between 00h jan 1st of current year and next week beginning
860         ELSE
861            idy = nday   + COUNT((/ sdjf%cltype == 'daily'                                 /))
862            imt = nmonth + COUNT((/ sdjf%cltype == 'monthly' .OR. idy > nmonth_len(nmonth) /))
863            iyr = nyear  + COUNT((/ sdjf%cltype == 'yearly'  .OR. imt == 13                /))
864            IF( idy > nmonth_len(nmonth) )   idy = 1
865            IF( imt == 13                )   imt = 1
866            isecwk = 0
867         ENDIF
868      ENDIF
869      !
870      ! find the last record to be read -> update sdjf%nreclast
871      indexyr = iyr - nyear + 1                 ! which  year are we looking for? previous(0), current(1) or next(2)?
872      indexmt = imt + 12 * ( indexyr - 1 )      ! which month are we looking for (relatively to current year)?
873      !
874      ! Last record to be read in the current file
875      ! Predefine the number of record in the file according of its type.
876      ! We could compare this number with the number of records in the file and make a stop if the 2 numbers do not match...
877      ! However this would be much less fexible (e.g. for tests) and will force to rewite input files according to nleapy...
878      IF    ( NINT(sdjf%freqh) == -12 ) THEN            ;   ireclast = 1    ! yearly mean: consider only 1 record
879      ELSEIF( NINT(sdjf%freqh) ==  -1 ) THEN                                ! monthly mean:
880         IF(     sdjf%cltype      == 'monthly' ) THEN   ;   ireclast = 1    !  consider that the file has  1 record
881         ELSE                                           ;   ireclast = 12   !  consider that the file has 12 record
882         ENDIF
883      ELSE                                                                  ! higher frequency mean (in hours)
884         IF(     sdjf%cltype      == 'monthly' ) THEN   ;   ireclast = NINT( 24. * REAL(nmonth_len(indexmt), wp) / sdjf%freqh )
885         ELSEIF( sdjf%cltype(1:4) == 'week'    ) THEN   ;   ireclast = NINT( 24. * 7.                            / sdjf%freqh )
886         ELSEIF( sdjf%cltype      == 'daily'   ) THEN   ;   ireclast = NINT( 24.                                 / sdjf%freqh )
887         ELSE                                           ;   ireclast = NINT( 24. * REAL( nyear_len(indexyr), wp) / sdjf%freqh )
888         ENDIF
889      ENDIF
890
891      sdjf%nreclast = ireclast
892      ! Allocate arrays for beginning/middle/end of each record (seconds since Jan. 1st 00h of nit000 year)
893      IF( ALLOCATED(sdjf%nrecsec) )   DEALLOCATE( sdjf%nrecsec )
894      ALLOCATE( sdjf%nrecsec( 0:ireclast ) )
895      !
896      IF    ( NINT(sdjf%freqh) == -12 ) THEN                                     ! yearly mean and yearly file
897         SELECT CASE( indexyr )
898         CASE(0)   ;   sdjf%nrecsec(0) = nsec1jan000 - nyear_len( 0 ) * idaysec
899         CASE(1)   ;   sdjf%nrecsec(0) = nsec1jan000
900         CASE(2)   ;   sdjf%nrecsec(0) = nsec1jan000 + nyear_len( 1 ) * idaysec
901         ENDSELECT
902         sdjf%nrecsec(1) = sdjf%nrecsec(0) + nyear_len( indexyr ) * idaysec
903      ELSEIF( NINT(sdjf%freqh) ==  -1 ) THEN                                     ! monthly mean:
904         IF(     sdjf%cltype      == 'monthly' ) THEN                            !    monthly file
905            sdjf%nrecsec(0   ) = nsec1jan000 + nmonth_beg(indexmt  )
906            sdjf%nrecsec(1   ) = nsec1jan000 + nmonth_beg(indexmt+1)
907         ELSE                                                                    !    yearly  file
908            ishift = 12 * ( indexyr - 1 )
909            sdjf%nrecsec(0:12) = nsec1jan000 + nmonth_beg(1+ishift:13+ishift)
910         ENDIF
911      ELSE                                                                       ! higher frequency mean (in hours)
912         IF(     sdjf%cltype      == 'monthly' ) THEN   ;   istart = nsec1jan000 + nmonth_beg(indexmt)
913         ELSEIF( sdjf%cltype(1:4) == 'week'    ) THEN   ;   istart = nsec1jan000 + isecwk
914         ELSEIF( sdjf%cltype      == 'daily'   ) THEN   ;   istart = nsec1jan000 + nmonth_beg(indexmt) + ( idy - 1 ) * idaysec
915         ELSEIF( indexyr          == 0         ) THEN   ;   istart = nsec1jan000 - nyear_len( 0 ) * idaysec
916         ELSEIF( indexyr          == 2         ) THEN   ;   istart = nsec1jan000 + nyear_len( 1 ) * idaysec
917         ELSE                                           ;   istart = nsec1jan000
918         ENDIF
919         zfreqs = sdjf%freqh * rhhmm * rmmss
920         DO jt = 0, sdjf%nreclast
921            sdjf%nrecsec(jt) = istart + NINT( zfreqs * REAL(jt,wp) )
922         END DO
923      ENDIF
924      !
925      IF( sdjf%ln_tint ) THEN   ! record time defined in the middle of the record, computed using an implementation
926                                ! of the rounded average that is valid over the full integer range
927         sdjf%nrecsec(1:sdjf%nreclast) = sdjf%nrecsec(0:sdjf%nreclast-1) / 2 + sdjf%nrecsec(1:sdjf%nreclast) / 2 + &
928            & MAX( MOD( sdjf%nrecsec(0:sdjf%nreclast-1), 2 ), MOD( sdjf%nrecsec(1:sdjf%nreclast), 2 ) )
929      END IF
930      !
931      sdjf%clname = fld_filename( sdjf, idy, imt, iyr )
932      !
933   END SUBROUTINE fld_def
934
935   
936   SUBROUTINE fld_clopn( sdjf )
937      !!---------------------------------------------------------------------
938      !!                    ***  ROUTINE fld_clopn  ***
939      !!
940      !! ** Purpose :   close/open the files
941      !!----------------------------------------------------------------------
942      TYPE(FLD)        , INTENT(inout) ::   sdjf       ! input field related variables
943      !
944      INTEGER, DIMENSION(2)  :: isave
945      LOGICAL  :: llprev, llnext, llstop
946      !!----------------------------------------------------------------------
947      !
948      llprev = sdjf%nrecsec(sdjf%nreclast) < nsec000_1jan000   ! file ends before the beginning of the job -> file may not exist
949      llnext = sdjf%nrecsec(       0     ) > nsecend_1jan000   ! file begins after the end of the job -> file may not exist
950
951      llstop = sdjf%ln_clim .OR. .NOT. ( llprev .OR. llnext )
952
953      IF( sdjf%num <= 0 .OR. .NOT. sdjf%ln_clim  ) THEN
954         IF( sdjf%num > 0 )   CALL iom_close( sdjf%num )   ! close file if already open
955         CALL iom_open( sdjf%clname, sdjf%num, ldstop = llstop, ldiof = LEN_TRIM(sdjf%wgtname) > 0 )
956      ENDIF
957      !
958      IF( sdjf%num <= 0 .AND. .NOT. llstop ) THEN   ! file not found but we do accept this...
959         !
960         IF( llprev ) THEN   ! previous file does not exist : go back to current and accept to read only the first record
961            CALL ctl_warn('previous file: '//TRIM(sdjf%clname)//' not found -> go back to current year/month/week/day file')
962            isave(1:2) = sdjf%nrecsec(sdjf%nreclast-1:sdjf%nreclast)   ! save previous file info
963            CALL fld_def( sdjf )   ! go back to current file
964            sdjf%nreclast = 1   ! force to use only the first record (do as if other were not existing...)
965            sdjf%nrecsec(0:1) = isave(1:2)
966         ENDIF
967         !
968         IF( llnext ) THEN   ! next     file does not exist : go back to current and accept to read only the last  record
969            CALL ctl_warn('next file: '//TRIM(sdjf%clname)//' not found -> go back to current year/month/week/day file')
970            isave(1:2) = sdjf%nrecsec(0:1)    ! save next file info
971            CALL fld_def( sdjf )   ! go back to current file
972            ! -> read last record but keep record info from the first record of next file
973            sdjf%nrecsec(sdjf%nreclast-1:sdjf%nreclast) = isave(1:2)
974            sdjf%nrecsec(0:sdjf%nreclast-2) = nflag
975         ENDIF
976         !
977         CALL iom_open( sdjf%clname, sdjf%num, ldiof = LEN_TRIM(sdjf%wgtname) > 0 )   
978         !
979      ENDIF
980      !
981   END SUBROUTINE fld_clopn
982
983
984   SUBROUTINE fld_fill( sdf, sdf_n, cdir, cdcaller, cdtitle, cdnam, knoprint )
985      !!---------------------------------------------------------------------
986      !!                    ***  ROUTINE fld_fill  ***
987      !!
988      !! ** Purpose :   fill the data structure (sdf) with the associated information
989      !!              read in namelist (sdf_n) and control print
990      !!----------------------------------------------------------------------
991      TYPE(FLD)  , DIMENSION(:)          , INTENT(inout) ::   sdf        ! structure of input fields (file informations, fields read)
992      TYPE(FLD_N), DIMENSION(:)          , INTENT(in   ) ::   sdf_n      ! array of namelist information structures
993      CHARACTER(len=*)                   , INTENT(in   ) ::   cdir       ! Root directory for location of flx files
994      CHARACTER(len=*)                   , INTENT(in   ) ::   cdcaller   ! name of the calling routine
995      CHARACTER(len=*)                   , INTENT(in   ) ::   cdtitle    ! description of the calling routine
996      CHARACTER(len=*)                   , INTENT(in   ) ::   cdnam      ! name of the namelist from which sdf_n comes
997      INTEGER                  , OPTIONAL, INTENT(in   ) ::   knoprint   ! no calling routine information printed
998      !
999      INTEGER  ::   jf   ! dummy indices
1000      !!---------------------------------------------------------------------
1001      !
1002      DO jf = 1, SIZE(sdf)
1003         sdf(jf)%clrootname = sdf_n(jf)%clname
1004         IF( TRIM(sdf_n(jf)%clname) /= 'NOT USED' )   sdf(jf)%clrootname = TRIM( cdir )//sdf(jf)%clrootname
1005         sdf(jf)%clname     = "not yet defined"
1006         sdf(jf)%freqh      = sdf_n(jf)%freqh
1007         sdf(jf)%clvar      = sdf_n(jf)%clvar
1008         sdf(jf)%ln_tint    = sdf_n(jf)%ln_tint
1009         sdf(jf)%ln_clim    = sdf_n(jf)%ln_clim
1010         sdf(jf)%cltype     = sdf_n(jf)%cltype
1011         sdf(jf)%num        = -1
1012         sdf(jf)%nbb        = 1  ! start with before data in 1
1013         sdf(jf)%naa        = 2  ! start with after  data in 2
1014         sdf(jf)%wgtname    = " "
1015         IF( LEN( TRIM(sdf_n(jf)%wname) ) > 0 )   sdf(jf)%wgtname = TRIM( cdir )//sdf_n(jf)%wname
1016         sdf(jf)%lsmname = " "
1017         IF( LEN( TRIM(sdf_n(jf)%lname) ) > 0 )   sdf(jf)%lsmname = TRIM( cdir )//sdf_n(jf)%lname
1018         sdf(jf)%vcomp      = sdf_n(jf)%vcomp
1019         sdf(jf)%rotn(:)    = .TRUE.   ! pretend to be rotated -> won't try to rotate data before the first call to fld_get
1020         IF( sdf(jf)%cltype(1:4) == 'week' .AND. nn_leapy == 0  )   &
1021            &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdf(jf)%clrootname)//') needs nn_leapy = 1')
1022         IF( sdf(jf)%cltype(1:4) == 'week' .AND. sdf(jf)%ln_clim )   &
1023            &   CALL ctl_stop('fld_clopn: weekly file ('//TRIM(sdf(jf)%clrootname)//') needs ln_clim = .FALSE.')
1024         sdf(jf)%nreclast   = -1 ! Set to non zero default value to avoid errors, is updated to meaningful value during fld_clopn
1025         sdf(jf)%igrd       = 0
1026         sdf(jf)%ibdy       = 0
1027         sdf(jf)%imap       => NULL()
1028         sdf(jf)%ltotvel    = .FALSE.
1029         sdf(jf)%lzint      = .FALSE.
1030      END DO
1031      !
1032      IF(lwp) THEN      ! control print
1033         WRITE(numout,*)
1034         IF( .NOT.PRESENT( knoprint) ) THEN
1035            WRITE(numout,*) TRIM( cdcaller )//' : '//TRIM( cdtitle )
1036            WRITE(numout,*) (/ ('~', jf = 1, LEN_TRIM( cdcaller ) ) /)
1037         ENDIF
1038         WRITE(numout,*) '   fld_fill : fill data structure with information from namelist '//TRIM( cdnam )
1039         WRITE(numout,*) '   ~~~~~~~~'
1040         WRITE(numout,*) '      list of files and frequency (>0: in hours ; <0 in months)'
1041         DO jf = 1, SIZE(sdf)
1042            WRITE(numout,*) '      root filename: '  , TRIM( sdf(jf)%clrootname ), '   variable name: ', TRIM( sdf(jf)%clvar )
1043            WRITE(numout,*) '         frequency: '      ,       sdf(jf)%freqh       ,   &
1044               &                  '   time interp: '    ,       sdf(jf)%ln_tint     ,   &
1045               &                  '   climatology: '    ,       sdf(jf)%ln_clim
1046            WRITE(numout,*) '         weights: '        , TRIM( sdf(jf)%wgtname    ),   &
1047               &                  '   pairing: '        , TRIM( sdf(jf)%vcomp      ),   &
1048               &                  '   data type: '      ,       sdf(jf)%cltype      ,   &
1049               &                  '   land/sea mask:'   , TRIM( sdf(jf)%lsmname    )
1050            call flush(numout)
1051         END DO
1052      ENDIF
1053      !
1054   END SUBROUTINE fld_fill
1055
1056
1057   SUBROUTINE wgt_list( sd, kwgt )
1058      !!---------------------------------------------------------------------
1059      !!                    ***  ROUTINE wgt_list  ***
1060      !!
1061      !! ** Purpose :   search array of WGTs and find a weights file entry,
1062      !!                or return a new one adding it to the end if new entry.
1063      !!                the weights data is read in and restructured (fld_weight)
1064      !!----------------------------------------------------------------------
1065      TYPE( FLD ), INTENT(in   ) ::   sd        ! field with name of weights file
1066      INTEGER    , INTENT(  out) ::   kwgt      ! index of weights
1067      !
1068      INTEGER ::   kw, nestid   ! local integer
1069      !!----------------------------------------------------------------------
1070      !
1071      !! search down linked list
1072      !! weights filename is either present or we hit the end of the list
1073      !
1074      !! because agrif nest part of filenames are now added in iom_open
1075      !! to distinguish between weights files on the different grids, need to track
1076      !! nest number explicitly
1077      nestid = 0
1078#if defined key_agrif
1079      nestid = Agrif_Fixed()
1080#endif
1081      DO kw = 1, nxt_wgt-1
1082         IF( ref_wgts(kw)%wgtname == sd%wgtname .AND. &
1083             ref_wgts(kw)%nestid  == nestid) THEN
1084            kwgt = kw
1085            RETURN
1086         ENDIF
1087      END DO
1088      kwgt = nxt_wgt
1089      CALL fld_weight( sd )
1090      !
1091   END SUBROUTINE wgt_list
1092
1093
1094   SUBROUTINE wgt_print( )
1095      !!---------------------------------------------------------------------
1096      !!                    ***  ROUTINE wgt_print  ***
1097      !!
1098      !! ** Purpose :   print the list of known weights
1099      !!----------------------------------------------------------------------
1100      INTEGER ::   kw   !
1101      !!----------------------------------------------------------------------
1102      !
1103      DO kw = 1, nxt_wgt-1
1104         WRITE(numout,*) 'weight file:  ',TRIM(ref_wgts(kw)%wgtname)
1105         WRITE(numout,*) '      ddims:  ',ref_wgts(kw)%ddims(1),ref_wgts(kw)%ddims(2)
1106         WRITE(numout,*) '     numwgt:  ',ref_wgts(kw)%numwgt
1107         WRITE(numout,*) '     jpiwgt:  ',ref_wgts(kw)%jpiwgt
1108         WRITE(numout,*) '     jpjwgt:  ',ref_wgts(kw)%jpjwgt
1109         WRITE(numout,*) '    botleft:  ',ref_wgts(kw)%botleft
1110         WRITE(numout,*) '   topright:  ',ref_wgts(kw)%topright
1111         IF( ref_wgts(kw)%cyclic ) THEN
1112            WRITE(numout,*) '       cyclical'
1113            IF( ref_wgts(kw)%overlap > 0 ) WRITE(numout,*) '              with overlap of ', ref_wgts(kw)%overlap
1114         ELSE
1115            WRITE(numout,*) '       not cyclical'
1116         ENDIF
1117         IF( ASSOCIATED(ref_wgts(kw)%data_wgt) )  WRITE(numout,*) '       allocated'
1118      END DO
1119      !
1120   END SUBROUTINE wgt_print
1121
1122
1123   SUBROUTINE fld_weight( sd )
1124      !!---------------------------------------------------------------------
1125      !!                    ***  ROUTINE fld_weight  ***
1126      !!
1127      !! ** Purpose :   create a new WGT structure and fill in data from file,
1128      !!              restructuring as required
1129      !!----------------------------------------------------------------------
1130      TYPE( FLD ), INTENT(in) ::   sd   ! field with name of weights file
1131      !!
1132      INTEGER ::   ji,jj,jn   ! dummy loop indices
1133      INTEGER ::   inum       ! local logical unit
1134      INTEGER ::   id         ! local variable id
1135      INTEGER ::   ipk        ! local vertical dimension
1136      INTEGER ::   zwrap      ! local integer
1137      LOGICAL ::   cyclical   !
1138      CHARACTER (len=5) ::   clname   !
1139      INTEGER , DIMENSION(4) ::   ddims
1140      INTEGER                ::   isrc
1141      REAL(wp), DIMENSION(jpi,jpj) ::   data_tmp
1142      !!----------------------------------------------------------------------
1143      !
1144      IF( nxt_wgt > tot_wgts ) THEN
1145        CALL ctl_stop("fld_weight: weights array size exceeded, increase tot_wgts")
1146      ENDIF
1147      !
1148      !! new weights file entry, add in extra information
1149      !! a weights file represents a 2D grid of a certain shape, so we assume that the current
1150      !! input data file is representative of all other files to be opened and processed with the
1151      !! current weights file
1152
1153      !! get data grid dimensions
1154      id = iom_varid( sd%num, sd%clvar, ddims )
1155
1156      !! now open the weights file
1157      CALL iom_open ( sd%wgtname, inum )   ! interpolation weights
1158      IF( inum > 0 ) THEN
1159
1160         !! determine whether we have an east-west cyclic grid
1161         !! from global attribute called "ew_wrap" in the weights file
1162         !! note that if not found, iom_getatt returns -999 and cyclic with no overlap is assumed
1163         !! since this is the most common forcing configuration
1164
1165         CALL iom_getatt(inum, 'ew_wrap', zwrap)
1166         IF( zwrap >= 0 ) THEN
1167            cyclical = .TRUE.
1168         ELSE IF( zwrap == -999 ) THEN
1169            cyclical = .TRUE.
1170            zwrap = 0
1171         ELSE
1172            cyclical = .FALSE.
1173         ENDIF
1174
1175         ref_wgts(nxt_wgt)%ddims(1) = ddims(1)
1176         ref_wgts(nxt_wgt)%ddims(2) = ddims(2)
1177         ref_wgts(nxt_wgt)%wgtname = sd%wgtname
1178         ref_wgts(nxt_wgt)%overlap = zwrap
1179         ref_wgts(nxt_wgt)%cyclic = cyclical
1180         ref_wgts(nxt_wgt)%nestid = 0
1181#if defined key_agrif
1182         ref_wgts(nxt_wgt)%nestid = Agrif_Fixed()
1183#endif
1184         !! weights file is stored as a set of weights (wgt01->wgt04 or wgt01->wgt16)
1185         !! for each weight wgtNN there is an integer array srcNN which gives the point in
1186         !! the input data grid which is to be multiplied by the weight
1187         !! they are both arrays on the model grid so the result of the multiplication is
1188         !! added into an output array on the model grid as a running sum
1189
1190         !! two possible cases: bilinear (4 weights) or bicubic (16 weights)
1191         id = iom_varid(inum, 'src05', ldstop=.FALSE.)
1192         IF( id <= 0 ) THEN   ;   ref_wgts(nxt_wgt)%numwgt = 4
1193         ELSE                 ;   ref_wgts(nxt_wgt)%numwgt = 16
1194         ENDIF
1195
1196         ALLOCATE( ref_wgts(nxt_wgt)%data_jpi(Nis0:Nie0,Njs0:Nje0,4) )
1197         ALLOCATE( ref_wgts(nxt_wgt)%data_jpj(Nis0:Nie0,Njs0:Nje0,4) )
1198         ALLOCATE( ref_wgts(nxt_wgt)%data_wgt(Nis0:Nie0,Njs0:Nje0,ref_wgts(nxt_wgt)%numwgt) )
1199
1200         DO jn = 1,4
1201            WRITE(clname,'(a3,i2.2)') 'src',jn
1202            CALL iom_get ( inum, jpdom_global, clname, data_tmp(:,:), cd_type = 'Z' )   !  no call to lbc_lnk
1203            DO_2D_00_00
1204!!$               isrc = NINT(data_tmp(ji,jj)) - 1
1205               isrc = INT(data_tmp(ji,jj)) - 1
1206               ref_wgts(nxt_wgt)%data_jpi(ji,jj,jn) = 1 + MOD(isrc,  ref_wgts(nxt_wgt)%ddims(1))
1207               ref_wgts(nxt_wgt)%data_jpj(ji,jj,jn) = 1 +     isrc / ref_wgts(nxt_wgt)%ddims(1)
1208            END_2D
1209         END DO
1210
1211         DO jn = 1, ref_wgts(nxt_wgt)%numwgt
1212            WRITE(clname,'(a3,i2.2)') 'wgt',jn
1213            CALL iom_get ( inum, jpdom_global, clname, data_tmp(:,:), cd_type = 'Z' )   !  no call to lbc_lnk
1214            DO_2D_00_00
1215               ref_wgts(nxt_wgt)%data_wgt(ji,jj,jn) = data_tmp(ji,jj)
1216            END_2D
1217         END DO
1218         CALL iom_close (inum)
1219 
1220         ! find min and max indices in grid
1221         ref_wgts(nxt_wgt)%botleft( 1) = MINVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1222         ref_wgts(nxt_wgt)%botleft( 2) = MINVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1223         ref_wgts(nxt_wgt)%topright(1) = MAXVAL(ref_wgts(nxt_wgt)%data_jpi(:,:,:))
1224         ref_wgts(nxt_wgt)%topright(2) = MAXVAL(ref_wgts(nxt_wgt)%data_jpj(:,:,:))
1225
1226         ! and therefore dimensions of the input box
1227         ref_wgts(nxt_wgt)%jpiwgt = ref_wgts(nxt_wgt)%topright(1) - ref_wgts(nxt_wgt)%botleft(1) + 1
1228         ref_wgts(nxt_wgt)%jpjwgt = ref_wgts(nxt_wgt)%topright(2) - ref_wgts(nxt_wgt)%botleft(2) + 1
1229
1230         ! shift indexing of source grid
1231         ref_wgts(nxt_wgt)%data_jpi(:,:,:) = ref_wgts(nxt_wgt)%data_jpi(:,:,:) - ref_wgts(nxt_wgt)%botleft(1) + 1
1232         ref_wgts(nxt_wgt)%data_jpj(:,:,:) = ref_wgts(nxt_wgt)%data_jpj(:,:,:) - ref_wgts(nxt_wgt)%botleft(2) + 1
1233
1234         ! create input grid, give it a halo to allow gradient calculations
1235         ! SA: +3 stencil is a patch to avoid out-of-bound computation in some configuration.
1236         ! a more robust solution will be given in next release
1237         ipk =  SIZE(sd%fnow, 3)
1238         ALLOCATE( ref_wgts(nxt_wgt)%fly_dta(ref_wgts(nxt_wgt)%jpiwgt+3, ref_wgts(nxt_wgt)%jpjwgt+3 ,ipk) )
1239         IF( ref_wgts(nxt_wgt)%cyclic ) ALLOCATE( ref_wgts(nxt_wgt)%col(1,ref_wgts(nxt_wgt)%jpjwgt+3,ipk) )
1240         !
1241         nxt_wgt = nxt_wgt + 1
1242         !
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, DIMENSION (ileni,ilenj), INTENT(in   ) :: zfieldn       ! array of forcing field with undeff for land points
1327      REAL, DIMENSION (ileni,ilenj), INTENT(  out) :: zfield        ! array of forcing field
1328      !
1329      REAL   , DIMENSION (ileni,ilenj)   :: zmat1, zmat2, zmat3, zmat4  ! local arrays
1330      REAL   , DIMENSION (ileni,ilenj)   :: zmat5, zmat6, zmat7, zmat8  !   -     -
1331      REAL   , DIMENSION (ileni,ilenj)   :: zlsm2d                      !   -     -
1332      REAL   , 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_00_00( 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_00_00( 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_00_00( 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_00_00( 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_00_00( 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%cltype /= '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%cltype /= 'yearly' )   WRITE(clname, '(a,"_m",i2.2)' ) TRIM( sdjf%clrootname ), kmonth   ! add month
1595      ENDIF
1596      IF(    sdjf%cltype == 'daily' .OR. sdjf%cltype(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%cltype(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.