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

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

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

Last change on this file since 1649 was 1628, checked in by smasson, 15 years ago

add daily type of forcing file, see ticket:377

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