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

source: branches/DEV_r1784_3DF/NEMO/OPA_SRC/SBC/fldread.F90 @ 1824

Last change on this file since 1824 was 1824, checked in by cbricaud, 14 years ago

correction for the case monthly file and fre=-1

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