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 branches/DEV_r1784_mid_year_merge_2010/NEMO/OPA_SRC/SBC – NEMO

source: branches/DEV_r1784_mid_year_merge_2010/NEMO/OPA_SRC/SBC/fldread.F90 @ 1970

Last change on this file since 1970 was 1970, checked in by acc, 14 years ago

ticket #684 step 5: Add in changes from the trunk between revisions 1821 and 1879.

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