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

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

source: trunk/NEMO/OPA_SRC/SBC/fldread.F90 @ 1275

Last change on this file since 1275 was 1275, checked in by rblod, 15 years ago

First introduction off interpolation off the fly, see ticket #279

  • Property svn:keywords set to Id
File size: 48.7 KB
Line 
1MODULE fldread
2   !!======================================================================
3   !!                       ***  MODULE  fldread  ***
4   !! Ocean forcing:  read input field for surface boundary condition
5   !!=====================================================================
6   !! History :  9.0  !  06-06  (G. Madec) Original code
7   !!                 !  05-08  (S. Alderson) Modified for Interpolation in memory
8   !!                 !         from input grid to model grid
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   fld_read      : read input fields used for the computation of the
13   !!                   surface boundary condition
14   !!----------------------------------------------------------------------
15   USE oce             ! ocean dynamics and tracers
16   USE dom_oce         ! ocean space and time domain
17   USE phycst          ! ???
18   USE daymod          ! calendar
19   USE in_out_manager  ! I/O manager
20   USE iom             ! I/O manager library
21   USE geo2ocean       ! for vector rotation on to model grid
22
23   IMPLICIT NONE
24   PRIVATE   
25
26   TYPE, PUBLIC ::   FLD_N      !: Namelist field informations
27      CHARACTER(len = 34) ::   clname      ! generic name of the NetCDF flux file
28      REAL(wp)            ::   freqh       ! frequency of each flux file
29      CHARACTER(len = 34) ::   clvar       ! generic name of the variable in the NetCDF flux file
30      LOGICAL             ::   ln_tint     ! time interpolation or not (T/F)
31      LOGICAL             ::   ln_clim     ! climatology or not (T/F)
32      CHARACTER(len = 7)  ::   cltype      ! type of data file 'monthly' or yearly'
33      CHARACTER(len = 34) ::   wname       ! generic name of a NetCDF weights file to be used, blank if not
34      CHARACTER(len = 34) ::   vcomp       ! symbolic component name if a vector that needs rotation
35                                           ! a string starting with "U" or "V" for each component   
36                                           ! chars 2 onwards identify which components go together 
37   END TYPE FLD_N
38
39   TYPE, PUBLIC ::   FLD        !: Input field related variables
40      CHARACTER(len = 256)            ::   clrootname   ! generic name of the NetCDF file
41      CHARACTER(len = 256)            ::   clname       ! current name of the NetCDF file
42      REAL(wp)                        ::   freqh        ! frequency of each flux file
43      CHARACTER(len = 34)             ::   clvar        ! generic name of the variable in the NetCDF flux file
44      LOGICAL                         ::   ln_tint      ! time interpolation or not (T/F)
45      LOGICAL                         ::   ln_clim      ! climatology or not (T/F)
46      CHARACTER(len = 7)              ::   cltype       ! type of data file 'monthly' or yearly'
47      INTEGER                         ::   num          ! iom id of the jpfld files to be read
48      REAL(wp)                        ::   swap_sec     ! swapping time in second since Jan. 1st 00h of nit000 year
49      REAL(wp) , DIMENSION(2)         ::   rec_b        ! before record (1: index, 2: second since Jan. 1st 00h of nit000 year
50      REAL(wp) , DIMENSION(2)         ::   rec_a        ! after  record (1: index, 2: second since Jan. 1st 00h of nit000 year
51      REAL(wp) , ALLOCATABLE, DIMENSION(:,:)   ::   fnow         ! input fields interpolated to now time step
52      REAL(wp) , ALLOCATABLE, DIMENSION(:,:,:) ::   fdta         ! 2 consecutive record of input fields
53      CHARACTER(len = 256)            ::   wgtname      ! current name of the NetCDF weight file acting as a key
54                                                        ! into the WGTLIST structure
55      CHARACTER(len = 34)             ::   vcomp        ! symbolic name for a vector component that needs rotation
56      LOGICAL ,  DIMENSION(2)         ::   rotn         ! flag to indicate whether field has been rotated
57   END TYPE FLD
58
59!$AGRIF_DO_NOT_TREAT
60
61   !! keep list of all weights variables so they're only read in once
62   !! need to add AGRIF directives not to process this structure
63   !! also need to force wgtname to include AGRIF nest number
64   TYPE         ::   WGT        !: Input weights related variables
65      CHARACTER(len = 256)                    ::   wgtname      ! current name of the NetCDF weight file
66      INTEGER , DIMENSION(2)                  ::   ddims        ! shape of input grid
67      INTEGER , DIMENSION(2)                  ::   botleft      ! top left corner of box in input grid containing
68                                                                ! current processor grid
69      INTEGER , DIMENSION(2)                  ::   topright     ! top right corner of box
70      INTEGER                                 ::   jpiwgt       ! width of box on input grid
71      INTEGER                                 ::   jpjwgt       ! height of box on input grid
72      INTEGER                                 ::   numwgt       ! number of weights (4=bilinear, 16=bicubic)
73      INTEGER                                 ::   nestid       ! for agrif, keep track of nest we're in
74      INTEGER                                 ::   offset       ! =0 when cyclic grid has coincident first/last columns,
75                                                                ! =1 when they assumed to be one grid spacing apart     
76                                                                ! =-1 otherwise
77      LOGICAL                                 ::   cyclic       ! east-west cyclic or not
78      INTEGER, DIMENSION(:,:,:), POINTER      ::   data_jpi     ! array of source integers
79      INTEGER, DIMENSION(:,:,:), POINTER      ::   data_jpj     ! array of source integers
80      REAL(wp), DIMENSION(:,:,:), POINTER     ::   data_wgt     ! array of weights on model grid
81      REAL(wp), DIMENSION(:,:), POINTER       ::   fly_dta      ! array of values on input grid
82      REAL(wp), DIMENSION(:,:), POINTER       ::   col2         ! temporary array for reading in columns
83   END TYPE WGT
84
85   INTEGER,     PARAMETER             ::   tot_wgts = 10
86   TYPE( WGT ), DIMENSION(tot_wgts)   ::   ref_wgts     ! array of wgts
87   INTEGER                            ::   nxt_wgt = 1  ! point to next available space in ref_wgts array
88
89!$AGRIF_END_DO_NOT_TREAT
90
91   PUBLIC   fld_read, fld_fill   ! called by sbc... modules
92
93   !!----------------------------------------------------------------------
94   !!   OPA 9.0 , LOCEAN-IPSL (2006)
95   !! $Id$
96   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
97   !!----------------------------------------------------------------------
98
99CONTAINS
100
101   SUBROUTINE fld_read( kt, kn_fsbc, sd )
102      !!---------------------------------------------------------------------
103      !!                    ***  ROUTINE fld_read  ***
104      !!                   
105      !! ** Purpose :   provide at each time step the surface ocean fluxes
106      !!                (momentum, heat, freshwater and runoff)
107      !!
108      !! ** Method  :   READ each input fields in NetCDF files using IOM
109      !!      and intepolate it to the model time-step.
110      !!         Several assumptions are made on the input file:
111      !!      blahblahblah....
112      !!----------------------------------------------------------------------
113      INTEGER  , INTENT(in   )               ::   kt        ! ocean time step
114      INTEGER  , INTENT(in   )               ::   kn_fsbc   ! sbc computation period (in time step)
115      TYPE(FLD), INTENT(inout), DIMENSION(:) ::   sd        ! input field related variables
116      !!
117      CHARACTER (LEN=34)                     ::   acomp     ! dummy weight name
118      INTEGER                                ::   kf, nf    ! dummy indices
119      INTEGER                                ::   imf       ! size of the structure sd
120      REAL(wp), DIMENSION(jpi,jpj)           ::   utmp, vtmp! temporary arrays for vector rotation
121
122      INTEGER  ::   jf         ! dummy indices
123      INTEGER  ::   kw         ! index into wgts array
124      REAL(wp) ::   zreclast   ! last record to be read in the current year file
125      REAL(wp) ::   zsecend    ! number of second since Jan. 1st 00h of nit000 year at nitend
126      LOGICAL  ::   llnxtyr    ! open next year file?
127      LOGICAL  ::   llstop     ! stop is the file is not existing
128      REAL(wp) ::   ztinta     ! ratio applied to after  records when doing time interpolation
129      REAL(wp) ::   ztintb     ! ratio applied to before records when doing time interpolation
130      CHARACTER(LEN=1000) ::   clfmt   ! write format
131      !!---------------------------------------------------------------------
132      !
133      imf = SIZE( sd )
134      !                                         ! ===================== !
135      DO jf = 1, imf                            !    LOOP OVER FIELD    !
136         !                                      ! ===================== !
137         !
138         IF( kt == nit000 )   CALL fld_init( sd(jf) )
139         !
140         ! read/update the after data?
141         IF( rsec_year + sec1jan000 > sd(jf)%swap_sec ) THEN
142
143            IF( sd(jf)%ln_tint ) THEN         ! time interpolation: swap before record field
144!CDIR COLLAPSE
145               sd(jf)%fdta(:,:,1) = sd(jf)%fdta(:,:,2)
146               sd(jf)%rotn(1)     = sd(jf)%rotn(2)
147            ENDIF
148
149            ! update record informations
150            CALL fld_rec( sd(jf) )
151
152            ! do we have to change the year/month of the forcing field??
153            IF( sd(jf)%ln_tint ) THEN
154               ! if we do time interpolation we will need to open next year/month file before the end of the current year/month
155               ! if it is the case, we are still before the end of the year/month when calling fld_rec so sd(jf)%rec_a(1) will
156               ! be larger than the record number that should be read for current year/month (for ex. 13 for monthly mean file)
157
158               ! last record to be read in the current file
159               IF( sd(jf)%freqh == -1. ) THEN             ;   zreclast = 12.
160               ELSE                             
161                  IF( sd(jf)%cltype == 'monthly' ) THEN   ;   zreclast = 24. / sd(jf)%freqh * REAL( nmonth_len(nmonth), wp )
162                  ELSE                                    ;   zreclast = 24. / sd(jf)%freqh * REAL( nyear_len(     1 ), wp )
163                  ENDIF
164               ENDIF
165             
166               ! do we need next year data?
167               IF( sd(jf)%rec_a(1) > zreclast ) THEN
168
169                  sd(jf)%rec_a(1) = 1.              ! force to read the first record of the next year
170
171                  IF( .NOT. sd(jf)%ln_clim ) THEN   ! close the current file and open a new one.
172                     
173                     llnxtyr  = sd(jf)%cltype /= 'monthly' .OR. nmonth == 12   ! do we need to open next year file?
174                     ! if the run finishes at the end of the current year/month, we do accept that next year/month file does
175                     ! not exist. If the run continue farther than the current year/month, next year/month file must exist
176                     zsecend = rsec_year + sec1jan000 + REAL(nitend - kt, wp) * rdttra(1)   ! second at the end of the run
177                     llstop = zsecend > sd(jf)%swap_sec                                 ! read more than 1 record of next year
178
179                     CALL fld_clopn( sd(jf), nyear + COUNT((/llnxtyr/)), nmonth + 1 - 12 * COUNT((/llnxtyr/)), llstop )
180
181                     IF( sd(jf)%num == 0 .AND. .NOT. llstop ) THEN    ! next year file is not existing
182                        CALL ctl_warn('next year/month file: '//TRIM(sd(jf)%clname)//' not existing -> back to current year/month')
183                        CALL fld_clopn( sd(jf), nyear, nmonth )       ! back to the current year/month
184                        sd(jf)%rec_a(1) = zreclast     ! force to read the last record to be read in the current year file
185                     ENDIF
186
187                  ENDIF
188               ENDIF
189       
190            ELSE
191               ! if we are not doing time interpolation, we must change the year/month of the file just afer switching
192               ! to the NEW year/month. If it is the case, we are at the beginning of the year/month when calling fld_rec
193               ! so sd(jf)%rec_a(1) = 1
194               IF( sd(jf)%rec_a(1) == 1 .AND. .NOT. sd(jf)%ln_clim )   CALL fld_clopn( sd(jf), nyear, nmonth )
195            ENDIF
196
197            ! read after data
198            IF( LEN(TRIM(sd(jf)%wgtname)) > 0 ) THEN
199               CALL wgt_list( sd(jf), kw )
200               CALL fld_interp( sd(jf)%num, sd(jf)%clvar, kw, sd(jf)%fdta(:,:,2), NINT( sd(jf)%rec_a(1) ) )
201            ELSE
202               CALL iom_get( sd(jf)%num, jpdom_data, sd(jf)%clvar, sd(jf)%fdta(:,:,2), NINT( sd(jf)%rec_a(1) ) )
203            ENDIF
204            sd(jf)%rotn(2) = .FALSE.
205
206         ENDIF
207         !                                      ! ===================== !
208      END DO                                    !  END LOOP OVER FIELD  !
209      !                                         ! ===================== !
210
211      IF( kt == nit000 .AND. lwp ) CALL wgt_print()
212
213      !! Vector fields may need to be rotated onto the local grid direction
214      !! This has to happen before the time interpolations
215      !! (sga: following code should be modified so that pairs arent searched for each time
216
217      DO jf = 1, imf
218         !! find vector rotations required
219         IF( LEN(TRIM(sd(jf)%vcomp)) > 0 ) THEN
220             !! east-west component has symbolic name starting with 'U'
221             IF( sd(jf)%vcomp(1:1) == 'U' ) THEN
222                !! found an east-west component, look for the north-south component
223                !! which has same symbolic name but with 'U' replaced with 'V'
224                nf = LEN_TRIM( sd(jf)%vcomp )
225                IF( nf == 1) THEN
226                   acomp = 'V'
227                ELSE
228                   acomp = 'V' // sd(jf)%vcomp(2:nf)
229                ENDIF
230                kf = -1
231                DO nf = 1, imf
232                  IF( TRIM(sd(nf)%vcomp) == TRIM(acomp) ) kf = nf
233                END DO
234                IF( kf > 0 ) THEN
235                   !! fields jf,kf are two components which need to be rotated together
236                   DO nf = 1,2
237                      !! check each time level of this pair
238                      IF( .NOT. sd(jf)%rotn(nf) .AND. .NOT. sd(kf)%rotn(nf) ) THEN
239                         utmp(:,:) = 0.0
240                         vtmp(:,:) = 0.0
241                         CALL rot_rep( sd(jf)%fdta(:,:,nf), sd(kf)%fdta(:,:,nf), 'T', 'en->i', utmp(:,:) )
242                         CALL rot_rep( sd(jf)%fdta(:,:,nf), sd(kf)%fdta(:,:,nf), 'T', 'en->j', vtmp(:,:) )
243                         sd(jf)%fdta(:,:,nf) = utmp(:,:)
244                         sd(kf)%fdta(:,:,nf) = vtmp(:,:)
245                         sd(jf)%rotn(nf) = .TRUE.
246                         sd(kf)%rotn(nf) = .TRUE.
247                         IF( lwp .AND. kt == nit000 ) &
248                                   WRITE(numout,*) 'fld_read: vector pair (',  &
249                                                   TRIM(sd(jf)%clvar),',',TRIM(sd(kf)%clvar), &
250                                                   ') rotated on to model grid'
251                      ENDIF
252                   END DO
253                ENDIF
254             ENDIF
255         ENDIF
256      END DO
257
258      !                                         ! ===================== !
259      DO jf = 1, imf                            !    LOOP OVER FIELD    !
260         !                                      ! ===================== !
261         !
262         ! update field at each kn_fsbc time-step
263         IF( MOD( kt-1, kn_fsbc ) == 0 ) THEN   
264            !
265            IF( sd(jf)%ln_tint ) THEN
266               IF(lwp .AND. kt - nit000 <= 100 ) THEN
267                  clfmt = "('fld_read: var ', a, ' kt = ', i8,' Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
268                     &    "' records b/a: ', i4.4, '/', i4.4, ' (', f7.2,'/', f7.2, ' days)')"
269                  WRITE(numout, clfmt)  TRIM( sd(jf)%clvar ), kt, nyear, nmonth, nday,   &
270                     & NINT(sd(jf)%rec_b(1)), NINT(sd(jf)%rec_a(1)), sd(jf)%rec_b(2)/rday, sd(jf)%rec_a(2)/rday
271               ENDIF
272               !
273               ztinta =  ( rsec_year + sec1jan000 - sd(jf)%rec_b(2) ) / ( sd(jf)%rec_a(2) - sd(jf)%rec_b(2) )
274               ztintb =  1. - ztinta
275!CDIR COLLAPSE
276               sd(jf)%fnow(:,:) = ztintb * sd(jf)%fdta(:,:,1) + ztinta * sd(jf)%fdta(:,:,2)
277            ELSE
278               IF(lwp .AND. kt - nit000 <= 100 ) THEN
279                  clfmt = "('fld_read: var ', a, ' kt = ', i8,' Y/M/D = ', i4.4,'/', i2.2,'/', i2.2," //   &
280                     &    "' record: ', i4.4, ' at ', f7.2, ' day')"
281                  WRITE(numout, clfmt)  TRIM(sd(jf)%clvar), kt, nyear, nmonth, nday, NINT(sd(jf)%rec_a(1)), sd(jf)%rec_a(2)/rday
282               ENDIF
283!CDIR COLLAPSE
284               sd(jf)%fnow(:,:) = sd(jf)%fdta(:,:,2)   ! piecewise constant field
285 
286            ENDIF
287            !
288         ENDIF
289
290         IF( kt == nitend )   CALL iom_close( sd(jf)%num )   ! Close the input files
291
292         !                                      ! ===================== !
293      END DO                                    !  END LOOP OVER FIELD  !
294      !                                         ! ===================== !
295   END SUBROUTINE fld_read
296
297
298   SUBROUTINE fld_init( sdjf )
299      !!---------------------------------------------------------------------
300      !!                    ***  ROUTINE fld_init  ***
301      !!
302      !! ** Purpose :  - if time interpolation, read before data
303      !!               - open current year file
304      !!
305      !! ** Method  :   
306      !!----------------------------------------------------------------------
307      TYPE(FLD), INTENT(inout) ::   sdjf        ! input field related variables
308      !!
309      LOGICAL :: llprevyr       ! are we reading previous year  file?
310      LOGICAL :: llprevmth      ! are we reading previous month file?
311      LOGICAL :: llprev         ! llprevyr .OR. llprevmth
312      INTEGER :: idvar          ! variable id
313      INTEGER :: inrec          ! number of record existing for this variable
314      INTEGER :: kwgt
315      CHARACTER(LEN=1000) ::   clfmt   ! write format
316      !!---------------------------------------------------------------------
317
318      ! some default definitions...
319      sdjf%num = 0   ! default definition for non-opened file
320      IF( sdjf%ln_clim )   sdjf%clname = TRIM( sdjf%clrootname )   ! file name defaut definition, never change in this case
321      llprevyr  = .FALSE.
322      llprevmth = .FALSE.
323           
324      ! define record informations
325      CALL fld_rec( sdjf )
326
327      IF( sdjf%ln_tint ) THEN ! we need to read the previous record and we will put it in the current record structure
328         
329         IF( sdjf%rec_b(1) == 0.e0 ) THEN   ! we redefine record sdjf%rec_b(1) with the last record of previous year file
330            IF( sdjf%freqh == -1. ) THEN    ! monthly mean
331               sdjf%rec_b(1) = 12.          ! force to read december mean
332            ELSE   
333               IF( sdjf%cltype == 'monthly' ) THEN   ! monthly file
334                  sdjf%rec_b(1) = 24. / sdjf%freqh * REAL( nmonth_len(nmonth-1), wp )   ! last record of previous month
335                  llprevmth = .NOT. sdjf%ln_clim                                        ! use previous month file?
336                  llprevyr  = .NOT. sdjf%ln_clim .AND. nmonth == 1                      ! use previous year  file?
337               ELSE                                  ! yearly file
338                  sdjf%rec_b(1) = 24. / sdjf%freqh * REAL( nyear_len(0), wp )           ! last record of year month
339                  llprevyr  = .NOT. sdjf%ln_clim                                        ! use previous year  file?
340               ENDIF
341            ENDIF
342         ENDIF
343         llprev = llprevyr .OR. llprevmth
344
345         CALL fld_clopn( sdjf, nyear - COUNT((/llprevyr/)), nmonth - COUNT((/llprevmth/)) + 12 * COUNT((/llprevyr/)), llprev )
346
347         ! if previous year/month file is not existing, we switch to the current year/month
348         IF( llprev .AND. sdjf%num == 0 ) THEN
349            CALL ctl_warn( 'previous year/month file: '//TRIM(sdjf%clname)//' not existing -> back to current year/month' )
350            ! we force to read the first record of the current year/month instead of last record of previous year/month
351            llprev = .false.
352            sdjf%rec_b(1) = 1.
353            CALL fld_clopn( sdjf, nyear, nmonth )
354         ENDIF
355         
356         IF( llprev ) THEN   ! check if the last record sdjf%rec_n(1) exists in the file
357            idvar = iom_varid( sdjf%num, sdjf%clvar )                                        ! id of the variable sdjf%clvar
358            IF( idvar <= 0 )   RETURN
359            inrec = iom_file( sdjf%num )%dimsz( iom_file( sdjf%num )%ndims(idvar), idvar )   ! size of the last dim of idvar
360            sdjf%rec_b(1) = MIN( sdjf%rec_b(1), REAL( inrec, wp ) )   ! make sure we select an existing record
361         ENDIF
362
363         ! read before data into sdjf%fdta(:,:,2) because we will swap data in the following part of fld_read
364         IF( LEN(TRIM(sdjf%wgtname)) > 0 ) THEN
365            CALL wgt_list( sdjf, kwgt )
366            CALL fld_interp( sdjf%num, sdjf%clvar, kwgt, sdjf%fdta(:,:,2), NINT( sdjf%rec_b(1) ) )
367         ELSE
368            CALL iom_get( sdjf%num, jpdom_data, sdjf%clvar, sdjf%fdta(:,:,2), NINT( sdjf%rec_b(1) ) )
369         ENDIF
370         sdjf%rotn(2) = .FALSE.
371
372         clfmt = "('fld_init : time-interpolation for ', a, ' read previous record = ', i4, ' at time = ', f7.2, ' days')"
373         IF(lwp) WRITE(numout, clfmt) TRIM(sdjf%clvar), NINT(sdjf%rec_b(1)), sdjf%rec_b(2)/rday
374
375         IF( llprev )   CALL iom_close( sdjf%num )   ! close previous year file (-> redefine sdjf%num to 0)
376
377      ENDIF
378
379      IF( sdjf%num == 0 )   CALL fld_clopn( sdjf, nyear, nmonth )   ! make sure current year/month file is opened
380
381      sdjf%swap_sec = rsec_year + sec1jan000 - 1.   ! force read/update the after data in the following part of fld_read
382     
383   END SUBROUTINE fld_init
384
385
386   SUBROUTINE fld_rec( sdjf )
387      !!---------------------------------------------------------------------
388      !!                    ***  ROUTINE fld_rec  ***
389      !!
390      !! ** Purpose :   compute rec_a, rec_b and swap_sec
391      !!
392      !! ** Method  :   
393      !!----------------------------------------------------------------------
394      TYPE(FLD), INTENT(inout) ::   sdjf        ! input field related variables
395      !!
396      INTEGER  ::   irec        ! record number
397      REAL(wp) ::   zrec        ! record number
398      REAL(wp) ::   ztmp        ! temporary variable
399      REAL(wp) ::   zfreq_sec   ! frequency mean (in seconds)
400      !!----------------------------------------------------------------------
401      !
402      IF( sdjf%freqh == -1. ) THEN      ! monthly mean
403         !
404         IF( sdjf%ln_tint ) THEN                 ! time interpolation, shift by 1/2 record
405            !
406            !                  INT( ztmp )
407            !                     /|\
408            !                    1 |    *----
409            !                    0 |----(             
410            !                      |----+----|--> time
411            !                      0   /|\   1   (nday/nmonth_len(nmonth))
412            !                           |   
413            !                           |   
414            !       forcing record :  nmonth
415            !                           
416            ztmp  = REAL( nday, wp ) / REAL( nmonth_len(nmonth), wp ) + 0.5
417         ELSE
418            ztmp  = 0.e0
419         ENDIF
420         irec = nmonth + INT( ztmp )
421
422         IF( sdjf%ln_tint ) THEN   ;   sdjf%swap_sec = rmonth_half(irec) + sec1jan000   ! swap at the middle of the month
423         ELSE                      ;   sdjf%swap_sec = rmonth_end( irec) + sec1jan000   ! swap at the end    of the month
424         ENDIF
425
426         sdjf%rec_a(:) = (/ REAL( irec, wp ), rmonth_half(irec) + sec1jan000 /)   ! define after  record number and time
427         irec = irec - 1                                                          ! move back to previous record
428         sdjf%rec_b(:) = (/ REAL( irec, wp ), rmonth_half(irec) + sec1jan000 /)   ! define before record number and time
429         !
430      ELSE                          ! higher frequency mean (in hours)
431         !
432         zfreq_sec = sdjf%freqh * 3600.   ! frequency mean (in seconds)
433         ! number of second since the beginning of the file
434         IF( sdjf%cltype == 'monthly' ) THEN   ;   ztmp = rsec_month   ! since Jan 1 of the current year
435         ELSE                                  ;   ztmp = rsec_year    ! since the first day of the current month
436         ENDIF
437         IF( sdjf%ln_tint ) THEN                ! time interpolation, shift by 1/2 record
438            !
439            !                  INT( ztmp )
440            !                     /|\
441            !                    2 |        *-----(
442            !                    1 |  *-----(
443            !                    0 |--(             
444            !                      |--+--|--+--|--+--|--> time
445            !                      0 /|\ 1 /|\ 2 /|\ 3 (rsec_year/zfreq_sec) or (rsec_month/zfreq_sec)
446            !                         |     |     |
447            !                         |     |     |
448            !       forcing record :  1     2     3
449            !                   
450            ztmp= ztmp / zfreq_sec + 0.5
451         ELSE                 
452            !
453            !                  INT( ztmp )
454            !                     /|\
455            !                    2 |           *-----(
456            !                    1 |     *-----(
457            !                    0 |-----(             
458            !                      |--+--|--+--|--+--|--> time
459            !                      0 /|\ 1 /|\ 2 /|\ 3 (rsec_year/zfreq_sec) or (rsec_month/zfreq_sec)
460            !                         |     |     |
461            !                         |     |     |
462            !       forcing record :  1     2     3
463            !                           
464            ztmp= ztmp / zfreq_sec
465         ENDIF
466         zrec = 1. + REAL( INT( ztmp ), wp )
467
468         ! after record index and second since Jan. 1st 00h of nit000 year
469         sdjf%rec_a(:) = (/ zrec, zfreq_sec * ( zrec - 0.5 ) + sec1jan000 /)
470         IF( sdjf%cltype == 'monthly' )   &   ! add the number of second between Jan 1 and the end of previous month
471            sdjf%rec_a(2) = sdjf%rec_a(2) + rday * REAL(SUM(nmonth_len(1:nmonth -1)), wp)   ! ok if nmonth=1
472
473         ! before record index and second since Jan. 1st 00h of nit000 year
474         zrec = zrec - 1.                           ! move back to previous record
475         sdjf%rec_b(:) = (/ zrec, zfreq_sec * ( zrec - 0.5 ) + sec1jan000 /)
476         IF( sdjf%cltype == 'monthly' )   &   ! add the number of second between Jan 1 and the end of previous month
477            sdjf%rec_b(2) = sdjf%rec_b(2) + rday * REAL(SUM(nmonth_len(1:nmonth -1)), wp)   ! ok if nmonth=1
478
479         ! swapping time in second since Jan. 1st 00h of nit000 year
480         IF( sdjf%ln_tint ) THEN   ;   sdjf%swap_sec =  sdjf%rec_a(2)                     ! swap at the middle of the record
481         ELSE                      ;   sdjf%swap_sec =  sdjf%rec_a(2) + 0.5 * zfreq_sec   ! swap at the end    of the record
482         ENDIF       
483         !
484      ENDIF
485      !
486   END SUBROUTINE fld_rec
487
488
489   SUBROUTINE fld_clopn( sdjf, kyear, kmonth, ldstop )
490      !!---------------------------------------------------------------------
491      !!                    ***  ROUTINE fld_clopn  ***
492      !!
493      !! ** Purpose :   update the file name and open the file
494      !!
495      !! ** Method  :   
496      !!----------------------------------------------------------------------
497      TYPE(FLD), INTENT(inout)           ::   sdjf     ! input field related variables
498      INTEGER  , INTENT(in   )           ::   kyear    ! year value
499      INTEGER  , INTENT(in   )           ::   kmonth   ! month value
500      LOGICAL  , INTENT(in   ), OPTIONAL ::   ldstop   ! stop if open to read a non-existing file (default = .TRUE.)
501
502      IF( sdjf%num /= 0 )   CALL iom_close( sdjf%num )   ! close file if already open
503      ! build the new filename if not climatological data
504      IF( .NOT. sdjf%ln_clim ) THEN   ;   WRITE(sdjf%clname, '(a,"_y",i4.4)' ) TRIM( sdjf%clrootname ), kyear    ! add year
505         IF( sdjf%cltype == 'monthly' )   WRITE(sdjf%clname, '(a,"m",i2.2)'  ) TRIM( sdjf%clname     ), kmonth   ! add month
506      ENDIF
507      CALL iom_open( sdjf%clname, sdjf%num, ldstop = ldstop )
508      !
509   END SUBROUTINE fld_clopn
510
511
512   SUBROUTINE fld_fill( sdf, sdf_n, cdir, cdcaller, cdtitle, cdnam )
513      !!---------------------------------------------------------------------
514      !!                    ***  ROUTINE fld_fill  ***
515      !!
516      !! ** Purpose :   fill sdf with sdf_n and control print
517      !!
518      !! ** Method  :   
519      !!----------------------------------------------------------------------
520      TYPE(FLD)  , DIMENSION(:), INTENT(inout) ::   sdf        ! structure of input fields (file informations, fields read)
521      TYPE(FLD_N), DIMENSION(:), INTENT(in   ) ::   sdf_n      ! array of namelist information structures
522      CHARACTER(len=*)         , INTENT(in   ) ::   cdir       ! Root directory for location of flx files
523      CHARACTER(len=*)         , INTENT(in   ) ::   cdcaller   !
524      CHARACTER(len=*)         , INTENT(in   ) ::   cdtitle    !
525      CHARACTER(len=*)         , INTENT(in   ) ::   cdnam      !
526      !
527      INTEGER  ::   jf       ! dummy indices
528      !!---------------------------------------------------------------------
529
530      DO jf = 1, SIZE(sdf)
531         sdf(jf)%clrootname = TRIM( cdir )//TRIM( sdf_n(jf)%clname )
532         sdf(jf)%freqh      = sdf_n(jf)%freqh
533         sdf(jf)%clvar      = sdf_n(jf)%clvar
534         sdf(jf)%ln_tint    = sdf_n(jf)%ln_tint
535         sdf(jf)%ln_clim    = sdf_n(jf)%ln_clim
536         IF( sdf(jf)%freqh == -1. ) THEN   ;    sdf(jf)%cltype = 'yearly'
537         ELSE                               ;    sdf(jf)%cltype = sdf_n(jf)%cltype
538         ENDIF
539         sdf(jf)%wgtname = " "
540         IF( LEN( TRIM(sdf_n(jf)%wname) ) > 0 )  &
541            sdf(jf)%wgtname = TRIM( cdir )//TRIM( sdf_n(jf)%wname )
542         sdf(jf)%vcomp   = sdf_n(jf)%vcomp
543      END DO
544
545      IF(lwp) THEN      ! control print
546         WRITE(numout,*)
547         WRITE(numout,*) TRIM( cdcaller )//' : '//TRIM( cdtitle )
548         WRITE(numout,*) (/ ('~', jf = 1, LEN_TRIM( cdcaller ) ) /)
549         WRITE(numout,*) '          '//TRIM( cdnam )//' Namelist'
550         WRITE(numout,*) '          list of files and frequency (>0: in hours ; <0 in months)'
551         DO jf = 1, SIZE(sdf)
552            WRITE(numout,*) '               root filename: '  , TRIM( sdf(jf)%clrootname ),   &
553               &                          ' variable name: '  , TRIM( sdf(jf)%clvar      )
554            WRITE(numout,*) '               frequency: '      ,       sdf(jf)%freqh       ,   &
555               &                          ' time interp: '    ,       sdf(jf)%ln_tint     ,   &
556               &                          ' climatology: '    ,       sdf(jf)%ln_clim     ,   &
557               &                          ' weights    : '    , TRIM( sdf(jf)%wgtname    ),   &
558               &                          ' pairing    : '    , TRIM( sdf(jf)%vcomp      ),   &
559               &                          ' data type: '      ,       sdf(jf)%cltype
560         END DO
561      ENDIF
562     
563   END SUBROUTINE fld_fill
564
565
566   SUBROUTINE wgt_list( sd, kwgt )
567      !!---------------------------------------------------------------------
568      !!                    ***  ROUTINE wgt_list  ***
569      !!
570      !! ** Purpose :   search array of WGTs and find a weights file
571      !!                entry, or return a new one adding it to the end
572      !!                if it is a new entry, the weights data is read in and
573      !!                restructured (fld_weight)
574      !!
575      !! ** Method  :   
576      !!----------------------------------------------------------------------
577      TYPE( FLD ),      INTENT(in)    ::   sd        ! field with name of weights file
578      INTEGER,      INTENT(inout)     ::   kwgt      ! index of weights
579      !!
580      INTEGER                         ::   kw
581      INTEGER                         ::   nestid
582      LOGICAL                         ::   found
583      !!----------------------------------------------------------------------
584      !
585      !! search down linked list
586      !! weights filename is either present or we hit the end of the list
587      found = .FALSE.
588
589      !! because agrif nest part of filenames are now added in iom_open
590      !! to distinguish between weights files on the different grids, need to track
591      !! nest number explicitly
592      nestid = 0
593#if defined key_agrif
594      nestid = Agrif_Fixed()
595#endif
596      DO kw = 1, nxt_wgt-1
597         IF( TRIM(ref_wgts(kw)%wgtname) == TRIM(sd%wgtname) .AND. &
598             ref_wgts(kw)%nestid == nestid) THEN
599            kwgt = kw
600            found = .TRUE.
601            EXIT
602         ENDIF
603      END DO
604      IF( .NOT.found ) THEN
605         kwgt = nxt_wgt
606         CALL fld_weight( sd )
607      ENDIF
608
609   END SUBROUTINE wgt_list
610
611   SUBROUTINE wgt_print( )
612      !!---------------------------------------------------------------------
613      !!                    ***  ROUTINE wgt_print  ***
614      !!
615      !! ** Purpose :   print the list of known weights
616      !!
617      !! ** Method  :   
618      !!----------------------------------------------------------------------
619      !!
620      INTEGER                         ::   kw
621      !!----------------------------------------------------------------------
622      !
623
624      DO kw = 1, nxt_wgt-1
625         WRITE(numout,*) 'weight file:  ',TRIM(ref_wgts(kw)%wgtname)
626         WRITE(numout,*) '      ddims:  ',ref_wgts(kw)%ddims(1),ref_wgts(kw)%ddims(2)
627         WRITE(numout,*) '     numwgt:  ',ref_wgts(kw)%numwgt
628         WRITE(numout,*) '     jpiwgt:  ',ref_wgts(kw)%jpiwgt
629         WRITE(numout,*) '     jpjwgt:  ',ref_wgts(kw)%jpjwgt
630         WRITE(numout,*) '    botleft:  ',ref_wgts(kw)%botleft
631         WRITE(numout,*) '   topright:  ',ref_wgts(kw)%topright
632         IF( ref_wgts(kw)%cyclic ) THEN
633            WRITE(numout,*) '       cyclical'
634            IF( ref_wgts(kw)%offset > 0 ) WRITE(numout,*) '                 with offset'
635         ELSE
636            WRITE(numout,*) '       not cyclical'
637         ENDIF
638         IF( ASSOCIATED(ref_wgts(kw)%data_wgt) )  WRITE(numout,*) '       allocated'
639      END DO
640
641   END SUBROUTINE wgt_print
642
643   SUBROUTINE fld_weight( sd )
644      !!---------------------------------------------------------------------
645      !!                    ***  ROUTINE fld_weight  ***
646      !!
647      !! ** Purpose :   create a new WGT structure and fill in data from 
648      !!                file, restructuring as required
649      !!
650      !! ** Method  :   
651      !!----------------------------------------------------------------------
652      TYPE( FLD ),      INTENT(in)            ::   sd            ! field with name of weights file
653      !!
654      INTEGER                                 ::   jn            ! dummy loop indices
655      INTEGER                                 ::   inum          ! temporary logical unit
656      INTEGER                                 ::   id            ! temporary variable id
657      CHARACTER (len=5)                       ::   aname
658      INTEGER , DIMENSION(3)                  ::   ddims
659      INTEGER , DIMENSION(jpi, jpj)           ::   data_src
660      REAL(wp), DIMENSION(jpi, jpj)           ::   data_tmp
661      REAL(wp), DIMENSION(:,:), ALLOCATABLE   ::   line2, lines  ! temporary array to read 2 lineumns
662      CHARACTER (len=34)                      ::   lonvar
663      LOGICAL                                 ::   cyclical
664      REAL(wp)                                ::   resid, dlon   ! temporary array to read 2 lineumns
665      INTEGER                                 ::   offset        ! temporary integer
666      !!----------------------------------------------------------------------
667      !
668      IF( nxt_wgt > tot_wgts ) THEN
669        CALL ctl_stop("fld_weights: weights array size exceeded, increase tot_wgts")
670      ENDIF
671      !
672      !! new weights file entry, add in extra information
673      !! a weights file represents a 2D grid of a certain shape, so we assume that the current
674      !! input data file is representative of all other files to be opened and processed with the
675      !! current weights file
676
677      !! open input data file (non-model grid)
678      CALL iom_open( sd%clname, inum )
679
680      !! get dimensions
681      id = iom_varid( inum, sd%clvar, ddims )
682
683      !! check for an east-west cyclic grid
684      !! try to guess name of longitude variable
685
686      lonvar = 'nav_lon'
687      id = iom_varid(inum, TRIM(lonvar), ldstop=.FALSE.)
688      IF( id <= 0 ) THEN
689         lonvar = 'lon'
690         id = iom_varid(inum, TRIM(lonvar), ldstop=.FALSE.)
691      ENDIF
692
693      offset = -1
694      cyclical = .FALSE.
695      IF( id > 0 ) THEN
696         !! found a longitude variable
697         !! now going to assume that grid is regular so we can read a single row
698
699         !! because input array is 2d, have to present iom with 2d array even though we only need 1d slice
700         !! worse, we cant pass line2(:,1) to iom_get since this is treated as a 1d array which doesnt match input file
701         ALLOCATE( lines(ddims(1),2) )
702         CALL iom_get(inum, jpdom_unknown, lonvar, lines(:,:), 1, kstart=(/1,1/), kcount=(/ddims(1),2/) )
703
704         !! find largest grid spacing
705         lines(1:ddims(1)-1,2) = lines(2:ddims(1),1) - lines(1:ddims(1)-1,1)
706         dlon = MAXVAL( lines(1:ddims(1)-1,2) )
707
708         resid = ABS(ABS(lines(ddims(1),1)-lines(1,1))-360.0)
709         IF( resid < rsmall ) THEN
710            !! end rows overlap in longitude
711            offset = 0
712            cyclical = .TRUE.
713         ELSEIF( resid < 2.0*dlon ) THEN
714            !! also call it cyclic if difference between end points is less than twice dlon from 360
715            offset = 1
716            cyclical = .TRUE.
717         ENDIF
718
719         DEALLOCATE( lines )
720
721      ELSE
722         !! guessing failed
723         !! read in first and last columns of data variable
724         !! since we dont know the name of the longitude variable (or even if there is one)
725         !! we assume that if these two columns are equal, file is cyclic east-west
726
727         !! because input array is 2d, have to present iom with 2d array even though we only need 1d slice
728         !! worse, we cant pass line2(1,:) to iom_get since this is treated as a 1d array which doesnt match input file
729         ALLOCATE( lines(2,ddims(2)), line2(2,ddims(2)) )
730         CALL iom_get(inum, jpdom_unknown, sd%clvar, line2(:,:), 1, kstart=(/1,1/), kcount=(/2,ddims(2)/) )
731         lines(2,:) = line2(1,:)
732
733         CALL iom_get(inum, jpdom_unknown, sd%clvar, line2(:,:), 1, kstart=(/ddims(1)-1,1/), kcount=(/2,ddims(2)/) )
734         lines(1,:) = line2(2,:)
735
736         resid = SUM( ABS(lines(1,:) - lines(2,:)) )
737         IF( resid < ddims(2)*rsmall ) THEN
738            offset = 0
739            cyclical = .TRUE.
740         ENDIF
741
742         DEALLOCATE( lines, line2 )
743      ENDIF
744
745      !! close it
746      CALL iom_close( inum )
747
748      !! now open the weights file
749
750      CALL iom_open ( sd%wgtname, inum )   ! interpolation weights
751      IF ( inum > 0 ) THEN
752
753         ref_wgts(nxt_wgt)%ddims(1) = ddims(1)
754         ref_wgts(nxt_wgt)%ddims(2) = ddims(2)
755         ref_wgts(nxt_wgt)%wgtname = sd%wgtname
756         ref_wgts(nxt_wgt)%offset = -1
757         ref_wgts(nxt_wgt)%cyclic = .FALSE.
758         IF( cyclical ) THEN
759            ref_wgts(nxt_wgt)%offset = offset
760            ref_wgts(nxt_wgt)%cyclic = .TRUE.
761         ENDIF
762         ref_wgts(nxt_wgt)%nestid = 0
763#if defined key_agrif
764         ref_wgts(nxt_wgt)%nestid = Agrif_Fixed()
765#endif
766         !! weights file is stored as a set of weights (wgt01->wgt04 or wgt01->wgt16)
767         !! for each weight wgtNN there is an integer array srcNN which gives the point in
768         !! the input data grid which is to be multiplied by the weight
769         !! they are both arrays on the model grid so the result of the multiplication is
770         !! added into an output array on the model grid as a running sum
771
772         !! two possible cases: bilinear (4 weights) or bicubic (16 weights)
773         id = iom_varid(inum, 'src05', ldstop=.FALSE.)
774         IF( id <= 0) THEN
775            ref_wgts(nxt_wgt)%numwgt = 4
776         ELSE
777            ref_wgts(nxt_wgt)%numwgt = 16
778         ENDIF
779
780         ALLOCATE( ref_wgts(nxt_wgt)%data_jpi(jpi,jpj,4) )
781         ALLOCATE( ref_wgts(nxt_wgt)%data_jpj(jpi,jpj,4) )
782         ALLOCATE( ref_wgts(nxt_wgt)%data_wgt(jpi,jpj,ref_wgts(nxt_wgt)%numwgt) )
783
784         DO jn = 1,4
785            aname = ' '
786            WRITE(aname,'(a3,i2.2)') 'src',jn
787            data_tmp(:,:) = 0
788            CALL iom_get ( inum, jpdom_unknown, aname, data_tmp(1:nlci,1:nlcj), &
789                           kstart=(/nimpp,njmpp/), kcount=(/nlci,nlcj/) )
790            data_src(:,:) = INT(data_tmp(:,:))
791            ref_wgts(nxt_wgt)%data_jpj(:,:,jn) = 1 + (data_src(:,:)-1) / ref_wgts(nxt_wgt)%ddims(1)
792            ref_wgts(nxt_wgt)%data_jpi(:,:,jn) = data_src(:,:) - ref_wgts(nxt_wgt)%ddims(1)*(ref_wgts(nxt_wgt)%data_jpj(:,:,jn)-1)
793         END DO
794
795         DO jn = 1, ref_wgts(nxt_wgt)%numwgt
796            aname = ' '
797            WRITE(aname,'(a3,i2.2)') 'wgt',jn
798            ref_wgts(nxt_wgt)%data_wgt(1:nlci,1:nlcj,jn) = 0.0
799            CALL iom_get ( inum, jpdom_unknown, aname, ref_wgts(nxt_wgt)%data_wgt(1:nlci,1:nlcj,jn), &
800                           kstart=(/nimpp,njmpp/), kcount=(/nlci,nlcj/) )
801         END DO
802         CALL iom_close (inum)
803 
804         ! find min and max indices in grid
805         ref_wgts(nxt_wgt)%botleft(1) = MINVAL(ref_wgts(nxt_wgt)%data_jpi(1:nlci,1:nlcj,:))
806         ref_wgts(nxt_wgt)%botleft(2) = MINVAL(ref_wgts(nxt_wgt)%data_jpj(1:nlci,1:nlcj,:))
807         ref_wgts(nxt_wgt)%topright(1) = MAXVAL(ref_wgts(nxt_wgt)%data_jpi(1:nlci,1:nlcj,:))
808         ref_wgts(nxt_wgt)%topright(2) = MAXVAL(ref_wgts(nxt_wgt)%data_jpj(1:nlci,1:nlcj,:))
809
810         ! and therefore dimensions of the input box
811         ref_wgts(nxt_wgt)%jpiwgt = ref_wgts(nxt_wgt)%topright(1) - ref_wgts(nxt_wgt)%botleft(1) + 1
812         ref_wgts(nxt_wgt)%jpjwgt = ref_wgts(nxt_wgt)%topright(2) - ref_wgts(nxt_wgt)%botleft(2) + 1
813
814         ! shift indexing of source grid
815         ref_wgts(nxt_wgt)%data_jpi(:,:,:) = ref_wgts(nxt_wgt)%data_jpi(:,:,:) - ref_wgts(nxt_wgt)%botleft(1) + 1
816         ref_wgts(nxt_wgt)%data_jpj(:,:,:) = ref_wgts(nxt_wgt)%data_jpj(:,:,:) - ref_wgts(nxt_wgt)%botleft(2) + 1
817
818         ! create input grid, give it a halo to allow gradient calculations
819         ALLOCATE( ref_wgts(nxt_wgt)%fly_dta(ref_wgts(nxt_wgt)%jpiwgt+2, ref_wgts(nxt_wgt)%jpjwgt+2) )
820         IF( ref_wgts(nxt_wgt)%cyclic ) ALLOCATE( ref_wgts(nxt_wgt)%col2(2,ref_wgts(nxt_wgt)%jpjwgt+2) )
821
822         nxt_wgt = nxt_wgt + 1
823
824      ELSE
825         CALL ctl_stop( '    fld_weight : unable to read the file ' )
826      ENDIF
827
828   END SUBROUTINE fld_weight
829
830   SUBROUTINE fld_interp(num, clvar, kw, dta, nrec)
831      !!---------------------------------------------------------------------
832      !!                    ***  ROUTINE fld_interp  ***
833      !!
834      !! ** Purpose :   apply weights to input gridded data to create data
835      !!                on model grid
836      !!
837      !! ** Method  :   
838      !!----------------------------------------------------------------------
839      INTEGER,          INTENT(in)                        ::   num                 ! stream number
840      CHARACTER(LEN=*), INTENT(in)                        ::   clvar               ! variable name
841      INTEGER,          INTENT(in)                        ::   kw                  ! weights number
842      REAL(wp),         INTENT(inout), DIMENSION(jpi,jpj) ::   dta                 ! output field on model grid
843      INTEGER,          INTENT(in)                        ::   nrec                ! record number to read (ie time slice)
844      !!
845      INTEGER, DIMENSION(2)                               ::   rec1,recn           ! temporary arrays for start and length
846      INTEGER                                             ::  jk, jn, jm           ! loop counters
847      INTEGER                                             ::  ni, nj               ! lengths
848      INTEGER                                             ::  jpimin,jpiwid        ! temporary indices
849      INTEGER                                             ::  jpjmin,jpjwid        ! temporary indices
850      INTEGER                                             ::  jpi1,jpi2,jpj1,jpj2  ! temporary indices
851      !!----------------------------------------------------------------------
852      !
853
854      !! for weighted interpolation we have weights at four corners of a box surrounding
855      !! a model grid point, each weight is multiplied by a grid value (bilinear case)
856      !! or by a grid value and gradients at the corner point (bicubic case)
857      !! so we need to have a 4 by 4 subgrid surrounding each model point to cover both cases
858
859      !! sub grid where we already have weights
860      jpimin = ref_wgts(kw)%botleft(1)
861      jpjmin = ref_wgts(kw)%botleft(2)
862      jpiwid = ref_wgts(kw)%jpiwgt
863      jpjwid = ref_wgts(kw)%jpjwgt
864
865      !! what we need to read into sub grid in order to calculate gradients
866      rec1(1) = MAX( jpimin-1, 1 )
867      rec1(2) = MAX( jpjmin-1, 1 )
868      recn(1) = MIN( jpiwid+2, ref_wgts(kw)%ddims(1)-rec1(1)+1 )
869      recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
870
871      !! where we need to read it to
872      jpi1 = 2 + rec1(1) - jpimin
873      jpj1 = 2 + rec1(2) - jpjmin
874      jpi2 = jpi1 + recn(1) - 1
875      jpj2 = jpj1 + recn(2) - 1
876
877      ref_wgts(kw)%fly_dta(:,:) = 0.0
878      CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%fly_dta(jpi1:jpi2,jpj1:jpj2), nrec, rec1, recn)
879
880      !! first four weights common to both bilinear and bicubic
881      !! note that we have to offset by 1 into fly_dta array because of halo
882      dta(:,:) = 0.0
883      DO jk = 1,4
884        DO jn = 1, jpj
885          DO jm = 1,jpi
886            ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
887            nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
888            dta(jm,jn) = dta(jm,jn) + ref_wgts(kw)%data_wgt(jm,jn,jk) * ref_wgts(kw)%fly_dta(ni+1,nj+1)
889          END DO
890        END DO
891      END DO
892
893      IF (ref_wgts(kw)%numwgt .EQ. 16) THEN
894
895        !! fix up halo points that we couldnt read from file
896        IF( jpi1 == 2 ) THEN
897           ref_wgts(kw)%fly_dta(jpi1-1,:) = ref_wgts(kw)%fly_dta(jpi1,:)
898        ENDIF
899        IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
900           ref_wgts(kw)%fly_dta(jpi2+1,:) = ref_wgts(kw)%fly_dta(jpi2,:)
901        ENDIF
902        IF( jpj1 == 2 ) THEN
903           ref_wgts(kw)%fly_dta(:,jpj1-1) = ref_wgts(kw)%fly_dta(:,jpj1)
904        ENDIF
905        IF( jpj2 + jpjmin - 1 == ref_wgts(kw)%ddims(2)+1 .AND. jpj2 .lt. jpjwid+2 ) THEN
906           ref_wgts(kw)%fly_dta(:,jpj2+1) = 2.0*ref_wgts(kw)%fly_dta(:,jpj2) - ref_wgts(kw)%fly_dta(:,jpj2-1)
907        ENDIF
908
909        !! if data grid is cyclic we can do better on east-west edges
910        !! but have to allow for whether first and last columns are coincident
911        IF( ref_wgts(kw)%cyclic ) THEN
912           rec1(2) = MAX( jpjmin-1, 1 )
913           recn(1) = 2
914           recn(2) = MIN( jpjwid+2, ref_wgts(kw)%ddims(2)-rec1(2)+1 )
915           jpj1 = 2 + rec1(2) - jpjmin
916           jpj2 = jpj1 + recn(2) - 1
917           IF( jpi1 == 2 ) THEN
918              rec1(1) = ref_wgts(kw)%ddims(1) - 1
919              CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col2(:,jpj1:jpj2), nrec, rec1, recn)
920              ref_wgts(kw)%fly_dta(jpi1-1,jpj1:jpj2) = ref_wgts(kw)%col2(ref_wgts(kw)%offset+1,jpj1:jpj2)
921           ENDIF
922           IF( jpi2 + jpimin - 1 == ref_wgts(kw)%ddims(1)+1 ) THEN
923              rec1(1) = 1
924              CALL iom_get( num, jpdom_unknown, clvar, ref_wgts(kw)%col2(:,jpj1:jpj2), nrec, rec1, recn)
925              ref_wgts(kw)%fly_dta(jpi2+1,jpj1:jpj2) = ref_wgts(kw)%col2(2-ref_wgts(kw)%offset,jpj1:jpj2)
926           ENDIF
927        ENDIF
928
929        ! gradient in the i direction
930        DO jk = 1,4
931          DO jn = 1, jpj
932            DO jm = 1,jpi
933              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
934              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
935              dta(jm,jn) = dta(jm,jn) + ref_wgts(kw)%data_wgt(jm,jn,jk+4) * 0.5 *         &
936                               (ref_wgts(kw)%fly_dta(ni+2,nj+1) - ref_wgts(kw)%fly_dta(ni,nj+1))
937            END DO
938          END DO
939        END DO
940
941        ! gradient in the j direction
942        DO jk = 1,4
943          DO jn = 1, jpj
944            DO jm = 1,jpi
945              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
946              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
947              dta(jm,jn) = dta(jm,jn) + ref_wgts(kw)%data_wgt(jm,jn,jk+8) * 0.5 *         &
948                               (ref_wgts(kw)%fly_dta(ni+1,nj+2) - ref_wgts(kw)%fly_dta(ni+1,nj))
949            END DO
950          END DO
951        END DO
952
953        ! gradient in the ij direction
954        DO jk = 1,4
955          DO jn = 1, jpj
956            DO jm = 1,jpi
957              ni = ref_wgts(kw)%data_jpi(jm,jn,jk)
958              nj = ref_wgts(kw)%data_jpj(jm,jn,jk)
959              dta(jm,jn) = dta(jm,jn) + ref_wgts(kw)%data_wgt(jm,jn,jk+12) * 0.25 * ( &
960                               (ref_wgts(kw)%fly_dta(ni+2,nj+2) - ref_wgts(kw)%fly_dta(ni  ,nj+2)) -   &
961                               (ref_wgts(kw)%fly_dta(ni+2,nj  ) - ref_wgts(kw)%fly_dta(ni  ,nj  )))
962            END DO
963          END DO
964        END DO
965
966      END IF
967
968   END SUBROUTINE fld_interp
969 
970END MODULE fldread
Note: See TracBrowser for help on using the repository browser.