source: branches/publications/ORCHIDEE-ICE_SurfaceMassBalance/src_driver/readdim2.f90 @ 8398

Last change on this file since 8398 was 4263, checked in by josefine.ghattas, 7 years ago
  • Added model and subroutines headings (comments) according to coding guidelines for stomate_io, dim2_driver, readdim2 and intersurf.
  • Added propset for src_paralllel
  • Property svn:keywords set to HeadURL Date Author Revision
File size: 104.1 KB
Line 
1!  ==============================================================================================================================\n
2!  MODULE       : readdim2
3!
4!  CONTACT      : orchidee-help _at_ ipsl.jussieu.fr
5!
6!  LICENCE      : IPSL (2006)
7!  This software is governed by the CeCILL licence see ORCHIDEE/ORCHIDEE_CeCILL.LIC
8!
9!>\BRIEF        This module contains subroutines for reading the forcing file for the dim2_driver.
10!!
11!!
12!!\n DESCRIPTION  : This module contains subroutines for reading the forcing file for the dim2_driver.
13!!                  Following subroutines are public and called from dim2_driver :
14!!                    - forcing_info : Open the forcing file and return information about the grid in the forcing file.
15!!                                     Prepare for a zoom if needed.
16!!                                     Initialization of parallelization related to the grid.
17!!                    - forcing_read : Return the forcing data for the current time step of the model. The forcing file will
18!!                                     be read if it has not already been done for the current time-step in the forcing file.
19!!                    - forcing_grid : Calculate the longitudes and latitudes of the model grid.
20!!
21!! RECENT CHANGE(S): None
22!!
23!! REFERENCE(S) : None
24!!   
25!! SVN     :
26!! $HeadURL$
27!! $Date$
28!! $Revision$
29!! \n
30!_ ================================================================================================================================
31
32MODULE readdim2
33
34   USE ioipsl_para
35   USE weather
36   USE TIMER
37   USE constantes
38   USE solar
39   USE grid
40   USE mod_orchidee_para
41
42   IMPLICIT NONE
43
44   PRIVATE
45   PUBLIC  :: forcing_read, forcing_info, forcing_grid
46
47   INTEGER, SAVE                            :: iim_full, jjm_full, llm_full, ttm_full
48   INTEGER, SAVE                            :: iim_zoom, jjm_zoom
49   INTEGER, SAVE                            :: iim_g_begin,jjm_g_begin,iim_g_end,jjm_g_end
50   REAL, SAVE, ALLOCATABLE, DIMENSION(:,:)  :: data_full, lon_full, lat_full
51   REAL, SAVE, ALLOCATABLE, DIMENSION(:)    :: lev_full
52   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: itau, i_index, j_index,j_index_g
53   INTEGER, SAVE                            :: i_test, j_test
54   INTEGER, SAVE                            :: printlev_loc   !! Local printlev
55   LOGICAL, SAVE                            :: allow_weathergen, interpol, daily_interpol
56   LOGICAL, SAVE, PUBLIC                    :: weathergen, is_watchout
57   REAL, SAVE                               :: merid_res, zonal_res
58   LOGICAL, SAVE                            :: have_zaxis=.FALSE.
59!-
60!- Heigh controls and data
61!-
62   LOGICAL, SAVE                            :: zfixed, zsigma, zhybrid, zlevels, zheight 
63   LOGICAL, SAVE                            :: zsamelev_uv 
64   REAL, SAVE                               :: zlev_fixed, zlevuv_fixed 
65   REAL, SAVE                               :: zhybrid_a, zhybrid_b 
66   REAL, SAVE                               :: zhybriduv_a, zhybriduv_b 
67
68CONTAINS
69
70!! ==============================================================================================================================\n
71!! SUBROUTINE   : forcing_info
72!!
73!>\BRIEF        Open the forcing file and return information about the grid in the forcing file.
74!!
75!!\n DESCRIPTION : This subroutine will get all the information from the forcing file and prepare for the zoom if needed.
76!!                 It returns to the caller the sizes of the data it will receive at the forcing_read call.
77!!                 This is important so that the caller can allocate the right space.
78!!
79!!
80!! RECENT CHANGE(S): None
81!!
82!! MAIN OUTPUT VARIABLE(S):
83!!
84!! REFERENCE(S) :
85!!
86!_ ================================================================================================================================
87  SUBROUTINE forcing_info(filename, iim, jjm, llm, tm, date0, dt_force, force_id)
88
89    IMPLICIT NONE
90
91!! 0.1 Input variables
92    CHARACTER(LEN=*), INTENT(in) :: filename     !! Name of the file to be opened
93
94!! 0.2 Output variables
95    INTEGER, INTENT(out)         :: iim          !! Size in x of the forcing data
96    INTEGER, INTENT(out)         :: jjm          !! Size in y of the forcing data
97    INTEGER, INTENT(out)         :: llm          !! Number of levels in the forcing data (not yet used)
98    INTEGER, INTENT(out)         :: tm           !! Time dimension of the forcing
99    REAL, INTENT(out)            :: date0        !! The date at which the forcing file starts (julian days)
100    REAL, INTENT(out)            :: dt_force     !! Time-step of the forcing file in seconds
101    INTEGER, INTENT(out)         :: force_id     !! Id of the forcing file
102
103!! 0.3 Local variables
104    CHARACTER(LEN=20)                  :: calendar_str
105    CHARACTER(LEN=200)                 :: printstr       !! temporary character string to contain error message
106    REAL                               :: juld_1, juld_2
107    REAL, ALLOCATABLE, DIMENSION(:,:)  :: fcontfrac
108    REAL, ALLOCATABLE, DIMENSION(:,:)  :: qair
109    LOGICAL                            :: contfrac_exists=.FALSE.
110    INTEGER                            :: NbPoint
111    INTEGER                            :: i_test,j_test
112    INTEGER                            :: i,j,ind,ttm_part
113    INTEGER, ALLOCATABLE, DIMENSION(:) :: index_l
114    REAL, ALLOCATABLE, DIMENSION(:,:)  :: lon, lat
115    REAL, ALLOCATABLE, DIMENSION(:)    :: lev, levuv
116
117    !-
118    CALL flininfo(filename,  iim_full, jjm_full, llm_full, ttm_full, force_id)
119    !-
120    IF ( printlev_loc>=3 ) WRITE(numout,*) 'forcing_info : Details from forcing file :', &
121         iim_full, jjm_full, llm_full, ttm_full
122    !-
123    IF ( llm_full < 1 ) THEN
124       have_zaxis = .FALSE.
125    ELSE
126       have_zaxis = .TRUE.
127    ENDIF
128    WRITE(numout,*) 'have_zaxis : ', llm_full, have_zaxis
129    !-
130    ttm_part = 2
131    ALLOCATE(itau(ttm_part))
132    ALLOCATE(data_full(iim_full, jjm_full),lon_full(iim_full, jjm_full),&
133         & lat_full(iim_full, jjm_full))
134    ALLOCATE(lev_full(llm_full))
135    ALLOCATE(fcontfrac(iim_full,jjm_full))
136    !-
137    lev_full(:) = zero
138    !-
139    dt_force=zero
140    CALL flinopen &
141         &  (filename, .FALSE., iim_full, jjm_full, llm_full, lon_full, lat_full, &
142         &   lev_full, ttm_part, itau, date0, dt_force, force_id)
143    IF ( dt_force == zero ) THEN
144       dt_force = itau(2) - itau(1) 
145       itau(:) = itau(:) / dt_force
146    ENDIF
147    !  WRITE(numout,*) "forcing_info : Forcing time step out of flinopen : ",dt_force
148    !-
149    !- What are the alowed options for the temportal interpolation
150    !-
151    !Config Key   = ALLOW_WEATHERGEN
152    !Config Desc  = Allow weather generator to create data
153    !Config If    = [-]
154    !Config Def   = n
155    !Config Help  = This flag allows the forcing-reader to generate
156    !Config         synthetic data if the data in the file is too sparse
157    !Config         and the temporal resolution would not be enough to
158    !Config         run the model.
159    !Config Units = [FLAG]
160    !-
161    allow_weathergen = .FALSE.
162    CALL getin_p('ALLOW_WEATHERGEN',allow_weathergen)
163    !-
164    !- The calendar was set by the forcing file. If no "calendar" attribute was
165    !- found then it is assumed to be gregorian,
166    !MM => FALSE !! it is NOT assumed anything !
167    !- else it is what ever is written in this attribute.
168    !-
169    CALL ioget_calendar(calendar_str)
170    i=INDEX(calendar_str,ACHAR(0))
171    IF ( i > 0 ) THEN
172       calendar_str(i:20)=' '
173    ENDIF
174    !  WRITE(numout,*) "forcing_info : Calendar used : ",calendar_str
175    IF ( calendar_str == 'XXXX' ) THEN
176       WRITE(numout,*) "forcing_info : The calendar was not found in the forcing file."
177       IF (allow_weathergen) THEN
178          ! Then change the calendar
179          CALL ioconf_calendar("noleap") 
180       ELSE
181          WRITE(numout,*) "forcing_info : We will force it to gregorian by default."
182          CALL ioconf_calendar("gregorian") !! = 365.2425 ; "noleap" = 365.0; "360d"; "julian"=365.25
183       ENDIF
184    ENDIF
185    WRITE(numout,*) "readdim2 : Calendar used by the model : ",calendar_str
186    IF (ttm_full .GE. 2) THEN
187       juld_1 = itau2date(itau(1), date0, dt_force)
188       juld_2 = itau2date(itau(2), date0, dt_force)
189    ELSE
190       juld_1 = 0
191       juld_2 = 0
192       CALL ipslerr_p ( 3, 'forcing_info','What is that only one time step in the forcing file ?', &
193            &         ' That can not be right.','verify forcing file.')
194    ENDIF
195    !-
196    !- Initialize one_year / one_day
197    CALL ioget_calendar (one_year, one_day)
198    !-
199    !- What is the distance between the two first states. From this we will deduce what is
200    !- to be done.
201    weathergen = .FALSE.
202    interpol = .FALSE.
203    daily_interpol = .FALSE.
204    is_watchout = .FALSE.
205    !-
206    IF ( ABS(ABS(juld_2-juld_1)-30.) .LE. 2.) THEN
207       IF ( allow_weathergen ) THEN
208          weathergen = .TRUE.
209          WRITE(numout,*) 'Using weather generator.' 
210       ELSE
211          CALL ipslerr_p ( 3, 'forcing_info', &
212               &         'This seems to be a monthly file.', &
213               &         'We should use a weather generator with this file.', &
214               &         'This should be allowed in the run.def')
215       ENDIF
216    ELSEIF (( ABS(juld_1-juld_2) .LE. 1./4.) .OR. ( ABS(juld_1-juld_2) .EQ. 1.)) THEN
217       interpol = .TRUE.
218       WRITE(numout,*) 'We will interpolate between the forcing data time steps.' 
219       IF ( ABS(juld_1-juld_2) .EQ. 1.) THEN
220          daily_interpol = .TRUE.
221       ENDIF
222    ELSE
223       ! Using the weather generator with data other than monthly ones probably
224       ! needs some thinking.
225       WRITE(numout,*) 'The time step is not suitable:',ABS(juld_1-juld_2),' days.'
226       CALL ipslerr_p ( 3, 'forcing_info','The time step is not suitable.', &
227            &         '','We cannot do anything with these forcing data.')
228    ENDIF
229    !-
230    !- redefine the forcing time step if the weather generator is activated
231    !-
232    IF ( weathergen ) THEN
233       !Config Key   = DT_WEATHGEN
234       !Config Desc  = Calling frequency of weather generator
235       !Config If    = ALLOW_WEATHERGEN
236       !Config Def   = 1800.
237       !Config Help  = Determines how often the weather generator
238       !Config         is called (time step in s). Should be equal
239       !Config         to or larger than Sechiba's time step (say,
240       !Config         up to 6 times Sechiba's time step or so).
241       !Config Units = [seconds]
242       dt_force = 1800.
243       CALL getin_p('DT_WEATHGEN',dt_force)
244    ENDIF
245    !-
246    !- Define the zoom
247    !-
248    !Config Key   = LIMIT_WEST
249    !Config Desc  = Western limit of region
250    !Config If    = [-]
251    !Config Def   = -180.
252    !Config Help  = Western limit of the region we are
253    !Config         interested in. Between -180 and +180 degrees
254    !Config         The model will use the smalest regions from
255    !Config         region specified here and the one of the forcing file.
256    !Config Units = [Degrees]
257    !-
258    limit_west = -180.
259    CALL getin_p('LIMIT_WEST',limit_west)
260    !-
261    !Config Key   = LIMIT_EAST
262    !Config Desc  = Eastern limit of region
263    !Config If    = [-]
264    !Config Def   = 180.
265    !Config Help  = Eastern limit of the region we are
266    !Config         interested in. Between -180 and +180 degrees
267    !Config         The model will use the smalest regions from
268    !Config         region specified here and the one of the forcing file.
269    !Config Units = [Degrees]
270    !-
271    limit_east = 180.
272    CALL getin_p('LIMIT_EAST',limit_east)
273    !-
274    !Config Key   = LIMIT_NORTH
275    !Config Desc  = Northern limit of region
276    !Config If    = [-]
277    !Config Def   = 90.
278    !Config Help  = Northern limit of the region we are
279    !Config         interested in. Between +90 and -90 degrees
280    !Config         The model will use the smalest regions from
281    !Config         region specified here and the one of the forcing file.
282    !Config Units = [Degrees]
283    !-
284    limit_north = 90.
285    CALL getin_p('LIMIT_NORTH',limit_north)
286    !-
287    !Config Key   = LIMIT_SOUTH
288    !Config Desc  = Southern limit of region
289    !Config If    = [-]
290    !Config Def   = -90.
291    !Config Help  = Southern limit of the region we are
292    !Config         interested in. Between 90 and -90 degrees
293    !Config         The model will use the smalest regions from
294    !Config         region specified here and the one of the forcing file.
295    !Config Units = [Degrees]
296    !-
297    limit_south = -90.
298    CALL getin_p('LIMIT_SOUTH',limit_south)
299    !-
300    !- Calculate domain size
301    !-
302    IF ( interpol ) THEN
303       !-
304       !- If we use temporal interpolation, then we cannot change the resolution (yet?)
305       !-
306       ALLOCATE(i_index(iim_full), j_index(jjm_full),j_index_g(jjm_full))
307       IF (is_root_prc) THEN
308
309          CALL domain_size (limit_west, limit_east, limit_north, limit_south,&
310               &         iim_full, jjm_full, lon_full, lat_full, iim_zoom, jjm_zoom,&
311               &         i_index, j_index_g)
312
313          j_index(:)=j_index_g(:)
314
315          ALLOCATE(qair(iim_full,jjm_full))
316          CALL flinget_buffer (force_id,'Qair',iim_full, jjm_full, 1, ttm_full,  1, 1, data_full)
317          CALL forcing_zoom(data_full, qair)
318
319          CALL flinquery_var(force_id, 'contfrac', contfrac_exists)
320          IF ( contfrac_exists ) THEN
321             WRITE(numout,*) "contfrac exist in the forcing file."
322             CALL flinget_buffer (force_id,'contfrac',iim_full, jjm_full, 1, ttm_full,  1, 1, data_full)
323             CALL forcing_zoom(data_full, fcontfrac)
324             WRITE(numout,*) "fcontfrac min/max :",MINVAL(fcontfrac(1:iim_zoom,1:jjm_zoom)),MAXVAL(fcontfrac(1:iim_zoom,1:jjm_zoom))
325          ELSE
326             fcontfrac(:,:)=1.
327          ENDIF
328
329
330          DO i=1,iim_zoom
331             DO j=1,jjm_zoom
332                IF ( fcontfrac(i,j) <= EPSILON(1.) ) THEN
333                   qair(i,j) = 999999.
334                ENDIF
335             ENDDO
336          ENDDO
337
338          ALLOCATE(index_l(iim_zoom*jjm_zoom))
339          !- In this point is returning the global NbPoint with the global index
340          CALL forcing_landind(iim_zoom,jjm_zoom,qair,.TRUE.,NbPoint,index_l,i_test,j_test)
341          !
342          ! Work out the vertical layers to be used
343          !
344          CALL forcing_vertical_ioipsl(force_id) 
345       ELSE
346          ALLOCATE(index_l(1))
347       ENDIF
348
349       ! Initiate global grid and parallelism
350       CALL bcast(iim_zoom)
351       CALL bcast(jjm_zoom)
352       CALL bcast(NbPoint)
353       CALL grid_set_glo(iim_zoom,jjm_zoom,NbPoint)
354       CALL grid_allocate_glo(4)
355
356       ! Check consistency in the number of mpi processors and the number of land points
357       ! in order to prevent an exception 
358       IF (NbPoint < mpi_size) THEN
359          WRITE(printstr,*) 'The number of landpoints found (', NbPoint, &
360               ') is less than the number of processors selected (', mpi_size,')' 
361          CALL ipslerr_p(3, 'forcing_info', 'Wrong parallelization options', &
362               TRIM(printstr), & 
363               'Decrease the number of processors for the current grid') 
364       ENDIF
365       
366       !
367       !- global index index_g is the index_l of root proc
368       IF (is_root_prc) index_g(:)=index_l(1:NbPoint)
369
370       DEALLOCATE(index_l)
371
372       !
373       ! Distribute to all processors the information on the forcing
374       !
375       CALL bcast(index_g)
376       CALL Init_orchidee_data_para_driver(nbp_glo,index_g)
377       CALL init_ioipsl_para
378
379       ! Initialize printlev_loc
380       printlev_loc=get_printlev('readdim2')
381       WRITE(numout,*) 'In readdim2, : standard PRINTLEV= ', printlev
382       WRITE(numout,*) 'In readdim2, : local PRINTLEV_readdim2= ', printlev_loc
383
384!     CALL Init_writeField_p
385
386       CALL bcast(jjm_zoom)
387       CALL bcast(i_index)
388       CALL bcast(j_index_g)
389       CALL bcast(zfixed) 
390       CALL bcast(zsigma) 
391       CALL bcast(zhybrid) 
392       CALL bcast(zlevels) 
393       CALL bcast(zheight) 
394       CALL bcast(zsamelev_uv) 
395       CALL bcast(zlev_fixed) 
396       CALL bcast(zlevuv_fixed) 
397       CALL bcast(zhybrid_a) 
398       CALL bcast(zhybrid_b) 
399       CALL bcast(zhybriduv_a) 
400       CALL bcast(zhybriduv_b) 
401       ind=0
402       DO j=1,jjm_zoom
403          IF ( (j >= jj_begin) .AND. (j <= jj_end) ) THEN
404             ind=ind+1
405             j_index(ind)=j_index_g(j)
406          ENDIF
407       ENDDO
408
409       jjm_zoom=jj_nb
410       iim_zoom=iim_g
411
412       !-
413       !- If we use the weather generator, then we read zonal and meridional resolutions
414       !- This should be unified one day...
415       !-
416    ELSEIF ( weathergen ) THEN
417       !-
418       !Config Key   = MERID_RES
419       !Config Desc  = North-South Resolution
420       !Config Def   = 2.
421       !Config If    = ALLOW_WEATHERGEN
422       !Config Help  = North-South Resolution of the region we are
423       !Config         interested in.
424       !Config Units = [Degrees]
425       !-
426       merid_res = 2.
427       CALL getin_p('MERID_RES',merid_res)
428       !-
429       !Config Key   = ZONAL_RES
430       !Config Desc  = East-West Resolution
431       !Config Def   = 2.
432       !Config If    = ALLOW_WEATHERGEN
433       !Config Help  = East-West Resolution of the region we are
434       !Config         interested in. In degrees
435       !Config Units = [Degrees]
436       !-
437       zonal_res = 2.
438       CALL getin_p('ZONAL_RES',zonal_res)
439       !-
440       !- Number of time steps is meaningless in this case
441       !-
442       !    ttm_full = HUGE( ttm_full )
443       !MM Number (realistic) of time steps for half hour dt
444       ttm_full = NINT(one_year * 86400. / dt_force)
445       !-
446       IF (is_root_prc) THEN
447
448          !- In this point is returning the global NbPoint with the global index
449          CALL weathgen_domain_size (limit_west,limit_east,limit_north,limit_south, &
450               zonal_res,merid_res,iim_zoom,jjm_zoom)
451          ALLOCATE(index_l(iim_zoom*jjm_zoom))
452       ENDIF
453       CALL bcast(iim_zoom)
454       CALL bcast(jjm_zoom)
455
456       ALLOCATE(lon(iim_zoom,jjm_zoom))
457       ALLOCATE(lat(iim_zoom,jjm_zoom))
458       ALLOCATE(lev(llm_full),levuv(llm_full))
459       
460       ! We need lon and lat now for weathgen_init
461       CALL forcing_grid (iim_zoom,jjm_zoom,llm_full,lon,lat,init_f=.TRUE.)
462       CALL forcing_vertical_ioipsl(-1) 
463
464       IF (is_root_prc) THEN
465          CALL weathgen_init &
466               &        (filename,dt_force,force_id,iim_zoom,jjm_zoom, &
467               &         zonal_res,merid_res,lon,lat,index_l,NbPoint)
468!!$,&
469!!$               &         i_index,j_index_g)
470       ELSE
471          ALLOCATE(index_l(1))
472          index_l(1)=1
473       ENDIF
474
475       CALL bcast(NbPoint)
476       CALL grid_set_glo(iim_zoom,jjm_zoom,NbPoint)
477       CALL grid_allocate_glo(4)
478
479       !
480       !- global index index_g is the index_l of root proc
481       IF (is_root_prc) index_g(:)=index_l(1:NbPoint)
482
483       DEALLOCATE(index_l)
484
485     CALL bcast(index_g)
486     CALL Init_orchidee_data_para_driver(nbp_glo,index_g)
487     CALL init_ioipsl_para
488!     CALL Init_writeField_p
489     CALL bcast(jjm_zoom)
490!!$       CALL bcast(i_index)
491!!$       CALL bcast(j_index_g)
492
493!!$       ind=0
494!!$       DO j=1,jjm_zoom
495!!$          IF ( (j >= jj_begin) .AND. (j <= jj_end) ) THEN
496!!$             ind=ind+1
497!!$             j_index(ind)=j_index_g(j)
498!!$          ENDIF
499!!$       ENDDO
500
501       jjm_zoom=jj_nb
502       iim_zoom=iim_g
503       !
504       CALL weathgen_read_file(force_id,iim_zoom,jjm_zoom)
505
506       !-
507    ELSE
508       !-
509       CALL ipslerr_p(3,'forcing_info','Neither interpolation nor weather generator is specified.','','')
510       !-
511    ENDIF
512    !-
513    !- Transfer the right information to the caller
514    !-
515    iim = iim_zoom
516    jjm = jjm_zoom
517    llm = 1
518    tm = ttm_full
519    !-
520    IF ( printlev_loc>=3 ) WRITE(numout,*) 'forcing_info : end : ', iim,jjm, llm,tm
521    !-
522  END SUBROUTINE forcing_info
523
524
525!! ==============================================================================================================================\n
526!! SUBROUTINE   : forcing_read
527!!
528!>\BRIEF        Return forcing data for the current time step
529!!
530!!\n DESCRIPTION : Return the forcing data for the current time step of the model. The forcing file will
531!!                 be read if it has not already been done for the current time-step in the forcing file.
532!!
533!! RECENT CHANGE(S): None
534!!
535!! MAIN OUTPUT VARIABLE(S):
536!!
537!! REFERENCE(S) :
538!!
539!_ ================================================================================================================================
540SUBROUTINE forcing_read &
541  & (filename, rest_id, lrstread, lrstwrite, &
542  &  itauin, istp, itau_split, split, nb_spread, netrad_cons, date0,   &
543  &  dt_force, iim, jjm, lon, lat, zlev, zlevuv, ttm, &
544  &  swdown, coszang, precip, snowf, tair, u, v, qair, pb, lwdown, &
545  &  fcontfrac, fneighbours, fresolution, &
546  &  SWnet, Eair, petAcoef, peqAcoef, petBcoef, peqBcoef, cdrag, ccanopy, &
547  &  kindex, nbindex, force_id)
548
549  IMPLICIT NONE
550
551    !! 0. Variable and parameter declaration
552    !! 0.1 Input variables
553   CHARACTER(LEN=*), INTENT(IN) :: filename   !! name of the file to be opened
554   INTEGER, INTENT(IN) :: force_id            !! FLINCOM file id. It is used to close the file at the end of the run.
555   INTEGER, INTENT(IN) :: rest_id             !! ID of restart file
556   LOGICAL, INTENT(IN) :: lrstread            !! read restart file?
557   LOGICAL, INTENT(IN) :: lrstwrite           !! write restart file?
558   INTEGER, INTENT(IN) :: itauin              !! time step for which we need the data
559   INTEGER, INTENT(IN) :: istp                !! time step for restart file
560   INTEGER, INTENT(IN) :: itau_split          !! Current step between 2 forcing times-step (it decides if it is time to read)
561   INTEGER, INTENT(IN) :: split               !! The number of time steps between two time-steps of the forcing
562   INTEGER, INTENT(IN) :: nb_spread           !! Over how many time steps do we spread the precipitation
563   LOGICAL, INTENT(IN) :: netrad_cons         !! flag that decides if net radiation should be conserved
564   REAL, INTENT(IN)    :: date0               !! The date at which the forcing file starts (julian days)
565   REAL, INTENT(IN)    :: dt_force            !! time-step of the forcing file in seconds
566   INTEGER, INTENT(IN) :: iim                 !! Size of the grid in x
567   INTEGER, INTENT(IN) :: jjm                 !! Size of the grid in y
568   INTEGER, INTENT(IN) :: ttm                 !! number of time steps in all in the forcing file
569   REAL,DIMENSION(iim,jjm), INTENT(IN) :: lon !! Longitudes
570   REAL,DIMENSION(iim,jjm), INTENT(IN) :: lat !! Latitudes
571
572    !! 0.2 Output variables
573   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: zlev        !! First Levels if it exists (ie if watchout file)
574   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: zlevuv      !! First Levels of the wind (equal precedent, if it exists)
575   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: swdown      !! Downward solar radiation (W/m^2)
576   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: coszang     !! Cosine of the solar zenith angle (unitless)
577   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: precip      !! Precipitation (Rainfall) (kg/m^2s)
578   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: snowf       !! Snowfall (kg/m^2s)
579   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: tair        !! 1st level (2m ? in off-line) air temperature (K)
580   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: u           !! 1st level (2m/10m ? in off-line) (in theory !) wind speed (m/s)
581   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: v           !! 1st level (2m/10m ? in off-line) (in theory !) wind speed (m/s)
582   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: qair        !! 1st level (2m ? in off-line) humidity (kg/kg)
583   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: pb          !! Surface pressure (Pa)
584   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: lwdown      !! Downward long wave radiation (W/m^2)
585   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: fcontfrac   !! Continental fraction (no unit)
586   REAL,DIMENSION(iim,jjm,2), INTENT(OUT) :: fresolution     !! resolution in x and y dimensions for each points
587   INTEGER,DIMENSION(iim,jjm,8), INTENT(OUT) :: fneighbours  !! land neighbours
588
589   !! From a WATCHOUT file :
590   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: SWnet       !! Net surface short-wave flux
591   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: Eair        !! Air potential energy
592   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: petAcoef    !! Coeficients A from the PBL resolution for T
593   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: peqAcoef    !! Coeficients A from the PBL resolution for q
594   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: petBcoef    !! Coeficients B from the PBL resolution for T
595   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: peqBcoef    !! Coeficients B from the PBL resolution for q
596   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: cdrag       !! Surface drag
597   REAL,DIMENSION(iim,jjm), INTENT(OUT) :: ccanopy     !! CO2 concentration in the canopy
598
599    !! 0.3 Modified variable
600   INTEGER, INTENT(INOUT) :: nbindex                   !! Number of land points
601   INTEGER,DIMENSION(iim*jjm), INTENT(INOUT) :: kindex !! Index of all land-points in the data (used for the gathering)
602
603    !! 0.4 Local variables
604   INTEGER :: ik,i,j
605
606   IF ( interpol ) THEN
607
608     CALL forcing_read_interpol &
609         (filename, itauin, itau_split, split, nb_spread, netrad_cons, date0,   &
610          dt_force, iim, jjm, lon, lat, zlev, zlevuv, ttm, &
611          swdown, coszang, precip, snowf, tair, u, v, qair, pb, lwdown, &
612          fcontfrac, fneighbours, fresolution, &
613          SWnet, Eair, petAcoef, peqAcoef, petBcoef, peqBcoef, cdrag, ccanopy, &
614          kindex, nbindex, force_id)
615
616   ELSEIF ( weathergen ) THEN
617
618      IF (lrstread) THEN
619         fcontfrac(:,:) = 1.0
620         WRITE(numout,*) 'contfrac : ', MINVAL(fcontfrac), MAXVAL(fcontfrac)
621      ENDIF
622
623      IF ( (itauin == 0).AND.(itau_split == 0) ) THEN
624         CALL weathgen_main (istp, istp, filename, force_id, iim, jjm, 1, &
625              rest_id, lrstread, lrstwrite, &
626              limit_west, limit_east, limit_north, limit_south, &
627              zonal_res, merid_res, lon, lat, date0, dt_force, &
628              kindex, nbindex, &
629              swdown, precip, snowf, tair, u, v, qair, pb, lwdown)
630      ELSE
631         CALL weathgen_main (itauin, istp, filename, force_id, iim, jjm, 1, &
632              rest_id, lrstread, lrstwrite, &
633              limit_west, limit_east, limit_north, limit_south, &
634              zonal_res, merid_res, lon, lat, date0, dt_force, &
635              kindex, nbindex, &
636              swdown, precip, snowf, tair, u, v, qair, pb, lwdown)
637      ENDIF
638
639      IF ( (itauin == 0).AND.(itau_split == 0) ) THEN
640         !---
641         !--- Allocate grid stuff
642         !---
643         CALL grid_init ( nbindex, 4, "RegLonLat", "ForcingGrid" )
644         !---
645         !--- Compute
646         !---
647         CALL grid_stuff(nbp_glo, iim_g, jjm_g, lon_g, lat_g, index_g)
648         !CALL grid_stuff (nbindex, iim, jjm, lon, lat, kindex)
649         DO ik=1,nbindex
650         
651            j = ((kindex(ik)-1)/iim) + 1
652            i = (kindex(ik) - (j-1)*iim)
653            !-
654            !- Store variable to help describe the grid
655            !- once the points are gathered.
656            !-
657            fneighbours(i,j,:) = neighbours(ik,:)
658            !
659            fresolution(i,j,:) = resolution(ik,:)
660         ENDDO
661      ENDIF
662   ELSE
663
664      CALL ipslerr_p(3,'forcing_read','Neither interpolation nor weather generator is specified.','','')
665
666   ENDIF
667
668   IF (.NOT. is_watchout) THEN
669      ! We have to compute Potential air energy
670      WHERE(tair(:,:) < val_exp) 
671         eair(:,:) = cp_air*tair(:,:)+cte_grav*zlev(:,:)
672      ENDWHERE
673   ENDIF
674
675END SUBROUTINE forcing_read
676
677!! ==============================================================================================================================\n
678!! SUBROUTINE   : forcing_read_interpol
679!!
680!>\BRIEF       
681!!
682!!\n DESCRIPTION :
683!!
684!! RECENT CHANGE(S): None
685!!
686!! MAIN OUTPUT VARIABLE(S):
687!!
688!! REFERENCE(S) :
689!!
690!_ ================================================================================================================================
691SUBROUTINE forcing_read_interpol &
692  & (filename, itauin, itau_split, split, nb_spread, netrad_cons, date0,   &
693  &  dt_force, iim, jjm, lon, lat, zlev, zlevuv, ttm, swdown, coszang, rainf, snowf, tair, &
694  &  u, v, qair, pb, lwdown, &
695  &  fcontfrac, fneighbours, fresolution, &
696  &  SWnet, Eair, petAcoef, peqAcoef, petBcoef, peqBcoef, cdrag, ccanopy, &
697  &  kindex, nbindex, force_id)
698!---------------------------------------------------------------------
699!- filename   : name of the file to be opened
700!- itauin     : time step for which we need the data
701!- itau_split : Where are we within the splitting
702!-              of the time-steps of the forcing files
703!-              (it decides IF we READ or not)
704!- split      : The number of time steps we do
705!-              between two time-steps of the forcing
706!- nb_spread  : Over how many time steps do we spread the precipitation
707!- netrad_cons: flag that decides if net radiation should be conserved.
708!- date0      : The date at which the forcing file starts (julian days)
709!- dt_force   : time-step of the forcing file in seconds
710!- iim        : Size of the grid in x
711!- jjm        : size of the grid in y
712!- lon        : Longitudes
713!- lat        : Latitudes
714!- zlev       : First Levels if it exists (ie if watchout file)
715!- zlevuv     : First Levels of the wind (equal precedent, if it exists)
716!- ttm        : number of time steps in all in the forcing file
717!- swdown     : Downward solar radiation (W/m^2)
718!- coszang    : Cosine of the solar zenith angle (unitless)
719!- rainf      : Rainfall (kg/m^2s)
720!- snowf      : Snowfall (kg/m^2s)
721!- tair       : 2m air temperature (K)
722!- u and v    : 2m (in theory !) wind speed (m/s)
723!- qair       : 2m humidity (kg/kg)
724!- pb         : Surface pressure (Pa)
725!- lwdown     : Downward long wave radiation (W/m^2)
726!- fcontfrac  : Continental fraction (no unit)
727!- fneighbours: land neighbours
728!- fresolution: resolution in x and y dimensions for each points
729!-
730!- From a WATCHOUT file :
731!- SWnet      : Net surface short-wave flux
732!- Eair       : Air potential energy
733!- petAcoef   : Coeficients A from the PBL resolution for T
734!- peqAcoef   : Coeficients A from the PBL resolution for q
735!- petBcoef   : Coeficients B from the PBL resolution for T
736!- peqBcoef   : Coeficients B from the PBL resolution for q
737!- cdrag      : Surface drag
738!- ccanopy    : CO2 concentration in the canopy
739!-
740!- kindex     : Index of all land-points in the data
741!-              (used for the gathering)
742!- nbindex    : Number of land points
743!- force_id   : FLINCOM file id.
744!-              It is used to close the file at the end of the run.
745!---------------------------------------------------------------------
746   IMPLICIT NONE
747!-
748   INTEGER,PARAMETER :: lm=1
749!-
750!- Input variables
751!-
752   CHARACTER(LEN=*) :: filename
753   INTEGER :: itauin, itau_split, split, nb_spread
754   LOGICAL :: netrad_cons
755   REAL    :: date0, dt_force
756   INTEGER :: iim, jjm, ttm
757   REAL,DIMENSION(:,:),INTENT(IN) :: lon, lat   !- LOCAL data array (size=iim,jjm)
758   INTEGER, INTENT(IN) :: force_id
759!-
760!- Output variables
761!-
762   REAL,DIMENSION(:,:),INTENT(OUT) :: zlev, zlevuv, &  !- LOCAL data array (size=iim,jjm)
763  &  swdown, coszang, rainf, snowf, tair, u, v, qair, pb, lwdown, &
764  &  fcontfrac
765   REAL,DIMENSION(:,:,:),INTENT(OUT) :: fresolution    !- LOCAL data array (size=iim,jjm,2)
766   INTEGER,DIMENSION(:,:,:),INTENT(OUT) :: fneighbours !- LOCAL data array (size=iim,jjm,8)
767   ! for watchout files
768   REAL,DIMENSION(:,:),INTENT(OUT) :: &
769  &  SWnet, Eair, petAcoef, peqAcoef, petBcoef, peqBcoef, cdrag, ccanopy
770   INTEGER,DIMENSION(:),INTENT(INOUT) :: kindex    !- LOCAL index of the map
771   INTEGER, INTENT(INOUT) :: nbindex
772!-
773!- Local variables
774!-
775   INTEGER, SAVE :: last_read=0
776   INTEGER, SAVE :: itau_read, itau_read_nm1=0, itau_read_n=0
777   REAL,SAVE,ALLOCATABLE,DIMENSION(:,:) :: &
778  &  zlev_nm1, zlevuv_nm1, swdown_nm1, rainf_nm1, snowf_nm1, tair_nm1, u_nm1, v_nm1, qair_nm1, & 
779  &  pb_nm1, lwdown_nm1, & 
780  &  zlev_n, zlevuv_n, swdown_n, rainf_n, snowf_n, tair_n, u_n, v_n, qair_n, &
781  &  pb_n, lwdown_n, mean_coszang
782
783   REAL,SAVE,ALLOCATABLE,DIMENSION(:,:) :: &
784  &  startday_n, startday_nm1, daylength_n, daylength_nm1, tmax_n, tmax_nm1, tmin_nm1, tmin_nm2, tmin_n,   &
785  &  qsatta, qsattmin_n, qsattmin_nm1, qmin_n, qmin_nm1, qmax_n, qmax_nm1, qsa
786   REAL,SAVE    :: hour
787
788!  just for grid stuff if the forcing file is a watchout file
789   REAL, ALLOCATABLE, DIMENSION(:,:) :: tmpdata
790   ! variables to be read in watchout files
791   REAL,SAVE,ALLOCATABLE,DIMENSION(:,:) :: &
792  &  SWnet_nm1, Eair_nm1, petAcoef_nm1, peqAcoef_nm1, petBcoef_nm1, peqBcoef_nm1, cdrag_nm1, ccanopy_nm1, &
793  &  SWnet_n, Eair_n, petAcoef_n, peqAcoef_n, petBcoef_n, peqBcoef_n, cdrag_n, ccanopy_n
794   REAL, SAVE :: julian_for ! Date of the forcing to be read
795   REAL :: julian, ss, rw
796!jur,
797   REAL, SAVE :: julian0    ! First day of this year
798   INTEGER :: yy, mm, dd, is, i, j, ik
799   REAL(r_std), DIMENSION(2) :: min_resol, max_resol
800!  if Wind_N and Wind_E are in the file (and not just Wind)
801   LOGICAL, SAVE :: wind_N_exists=.FALSE.
802   LOGICAL       :: wind_E_exists=.FALSE.
803   LOGICAL, SAVE :: contfrac_exists=.FALSE.
804   LOGICAL, SAVE :: neighbours_exists=.FALSE.
805   LOGICAL, SAVE :: initialized = .FALSE.
806   LOGICAL :: check=.FALSE.
807!  to bypass FPE problem on integer convertion with missing_value heigher than precision
808   INTEGER, PARAMETER                         :: undef_int = 999999999
809!---------------------------------------------------------------------
810
811   itau_read = MOD((itauin-1),ttm)+1
812
813   IF (check) THEN
814      WRITE(numout,*) &
815           " FORCING READ : itauin, itau_read, itau_split : ",&
816           itauin, itau_read, itau_split
817   ENDIF
818
819!-
820!- This part initializes the reading of the forcing. As you can see
821!- we only go through here if both time steps are zero.
822!-
823   IF ( (itau_read == 0).AND.(itau_split == 0) ) THEN
824!-
825!- Tests on forcing file type
826     CALL flinquery_var(force_id, 'Wind_N', wind_N_exists)
827     IF ( wind_N_exists ) THEN
828        CALL flinquery_var(force_id, 'Wind_E', wind_E_exists)
829        IF ( .NOT. wind_E_exists ) THEN
830           CALL ipslerr_p(3,'forcing_read_interpol', &
831   &             'Variable Wind_E does not exist in forcing file', &
832   &             'But variable Wind_N exists.','Please, rename Wind_N to Wind;') 
833        ENDIF
834     ENDIF
835     CALL flinquery_var(force_id, 'levels', is_watchout)
836     IF ( is_watchout ) THEN
837        WRITE(numout,*) "Read a Watchout File."
838     ENDIF
839     CALL flinquery_var(force_id, 'contfrac', contfrac_exists)
840!-
841     IF (check) WRITE(numout,*) 'ALLOCATE all the memory needed'
842!-
843     ALLOCATE &
844    &  (swdown_nm1(iim,jjm), rainf_nm1(iim,jjm), snowf_nm1(iim,jjm), &
845    &   tair_nm1(iim,jjm), u_nm1(iim,jjm), v_nm1(iim,jjm), qair_nm1(iim,jjm), &
846    &   pb_nm1(iim,jjm), lwdown_nm1(iim,jjm))
847     ALLOCATE &
848    &  (swdown_n(iim,jjm), rainf_n(iim,jjm), snowf_n(iim,jjm), &
849    &   tair_n(iim,jjm), u_n(iim,jjm), v_n(iim,jjm), qair_n(iim,jjm), &
850    &   pb_n(iim,jjm), lwdown_n(iim,jjm))
851
852     IF(daily_interpol) THEN
853        ALLOCATE &
854       &  (startday_n(iim,jjm), startday_nm1(iim,jjm), daylength_n(iim,jjm), &
855       &   daylength_nm1(iim,jjm), tmax_n(iim,jjm), tmax_nm1(iim,jjm), tmin_n(iim,jjm), &
856       &   tmin_nm1(iim,jjm), tmin_nm2(iim,jjm), qsatta(iim,jjm), qsattmin_n(iim,jjm), qsattmin_nm1(iim,jjm), &
857       &   qmin_n(iim,jjm), qmin_nm1(iim,jjm), qmax_n(iim,jjm), qmax_nm1(iim,jjm), qsa(iim,jjm) )
858     ENDIF
859
860
861     ALLOCATE &
862    &  (zlev_nm1(iim,jjm), zlev_n(iim,jjm), zlevuv_nm1(iim,jjm), zlevuv_n(iim,jjm), &
863    &   SWnet_nm1(iim,jjm), Eair_nm1(iim,jjm), cdrag_nm1(iim,jjm), ccanopy_nm1(iim,jjm), &
864    &   petAcoef_nm1(iim,jjm), peqAcoef_nm1(iim,jjm), petBcoef_nm1(iim,jjm), peqBcoef_nm1(iim,jjm), &
865    &   SWnet_n(iim,jjm), Eair_n(iim,jjm), cdrag_n(iim,jjm), ccanopy_n(iim,jjm), &
866    &   petAcoef_n(iim,jjm), peqAcoef_n(iim,jjm), petBcoef_n(iim,jjm), peqBcoef_n(iim,jjm))
867     ALLOCATE &
868    &  (mean_coszang(iim,jjm))
869!-
870     IF (check) WRITE(numout,*) 'Memory ALLOCATED'
871!-
872!- We need for the driver surface air temperature and humidity before the
873!- the loop starts. Thus we read it here.
874!-     
875     CALL forcing_just_read (iim, jjm, zlev, zlevuv, ttm, 1, 1, &
876          &  swdown, rainf, snowf, tair, &
877          &  u, v, qair, pb, lwdown, &
878          &  SWnet, Eair, petAcoef, peqAcoef, petBcoef, peqBcoef, cdrag, ccanopy, &
879          &  force_id, wind_N_exists, check)
880!----
881
882!-- First in case it's not a watchout file
883     IF ( .NOT. is_watchout ) THEN
884        IF ( contfrac_exists ) THEN
885           WRITE(numout,*) "contfrac exist in the forcing file."
886           CALL flinget_buffer (force_id,'contfrac',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
887           CALL forcing_zoom(data_full, fcontfrac)
888           WRITE(numout,*) "fcontfrac min/max :",MINVAL(fcontfrac(1:iim_zoom,jjm_zoom)),MAXVAL(fcontfrac(1:iim_zoom,jjm_zoom))
889           !
890           ! We need to make sure that when we gather the points we pick all
891           ! the points where contfrac is above 0. Thus we prepare tair for
892           ! subroutine forcing_landind
893           !
894           DO i=1,iim
895              DO j=1,jjm
896                 IF ( j==1 .AND. i<ii_begin) fcontfrac(i,j)=0.     ! bande de recouvrement du scatter2D
897                 IF ( j==jjm .AND. i>ii_end) fcontfrac(i,j)=0.     ! => on mets les données qu'on veut pas du noeud à missing_value
898                 IF ( fcontfrac(i,j) <= EPSILON(1.) ) THEN
899                    tair(i,j) = 999999.
900                 ENDIF
901              ENDDO
902           ENDDO
903        ELSE
904           fcontfrac(:,:) = 1.0
905        ENDIF
906        !---
907        !--- Create the index table
908        !---
909        !--- This job return a LOCAL kindex
910        CALL forcing_landind(iim, jjm, tair, check, nbindex, kindex, i_test, j_test)
911#ifdef CPP_PARA
912        ! We keep previous function forcing_landind, only to get a valid (i_test,j_test)
913        ! Force nbindex points for parallel computation
914        nbindex=nbp_loc
915        CALL scatter(index_g,kindex(1:nbindex))
916        kindex(1:nbindex)=kindex(1:nbindex)-(jj_begin-1)*iim_g
917#endif
918        ik=MAX(nbindex/2,1)
919        j_test = (((kindex(ik)-1)/iim) + 1)
920        i_test = (kindex(ik) - (j_test-1)*iim)
921        IF (check) THEN
922           WRITE(numout,*) 'New test point is : ', i_test, j_test
923        ENDIF
924        !---
925        !--- Allocate grid stuff
926        !---
927        CALL grid_init ( nbindex, 4, "RegLonLat", "ForcingGrid" )
928        !---
929        !--- All grid variables
930        !---
931        CALL grid_stuff(nbp_glo, iim_g, jjm_g, lon_g, lat_g, index_g)
932        DO ik=1,nbindex
933           !
934           j = ((kindex(ik)-1)/iim) + 1
935           i = (kindex(ik) - (j-1)*iim)
936           !-
937           !- Store variable to help describe the grid
938           !- once the points are gathered.
939              !-
940           fneighbours(i,j,:) = neighbours(ik,:)
941           !
942           fresolution(i,j,:) = resolution(ik,:)
943        ENDDO
944     ELSE
945!-- Second, in case it is a watchout file
946        ALLOCATE (tmpdata(iim,jjm))
947        tmpdata(:,:) = 0.0
948!--
949        IF ( .NOT. contfrac_exists ) THEN
950           CALL ipslerr_p (3,'forcing_read_interpol', &
951 &          'Could get contfrac variable in a watchout file :',TRIM(filename), &
952 &          '(Problem with file ?)')
953        ENDIF
954        CALL flinget_buffer (force_id,'contfrac',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
955        CALL forcing_zoom(data_full, fcontfrac)
956        !
957        ! We need to make sure that when we gather the points we pick all
958        ! the points where contfrac is above 0. Thus we prepare tair for
959        ! subroutine forcing_landind
960        !
961        DO i=1,iim
962           DO j=1,jjm
963              IF ( j==1 .AND. i<ii_begin) fcontfrac(i,j)=0.
964              IF ( j==jjm .AND. i>ii_end) fcontfrac(i,j)=0.
965              IF ( fcontfrac(i,j) <= EPSILON(1.) ) THEN
966                 tair(i,j) = 999999.
967              ENDIF
968           ENDDO
969        ENDDO
970        !---
971        !--- Create the index table
972        !---
973        !--- This job return a LOCAL kindex
974        CALL forcing_landind(iim, jjm, tair, check, nbindex, kindex, i_test, j_test)
975#ifdef CPP_PARA
976        ! We keep previous function forcing_landind, only to get a valid (i_test,j_test)
977        ! Force nbindex points for parallel computation
978        nbindex=nbp_loc
979        CALL scatter(index_g,kindex)
980        kindex(:)=kindex(:)-offset
981!        kindex(:)=kindex(:)-(jj_begin-1)*iim_g
982#endif
983        ik=MAX(nbindex/2,1)
984        j_test = (((kindex(ik)-1)/iim) + 1)
985        i_test = (kindex(ik) - (j_test-1)*iim)
986        IF (check) THEN
987           WRITE(numout,*) 'New test point is : ', i_test, j_test
988        ENDIF
989        !---
990        !--- Allocate grid stuff
991        !---
992        CALL grid_init ( nbindex, 4, "RegLonLat", "ForcingGrid" )
993        neighbours(:,:) = -1
994        resolution(:,:) = 0.
995        min_resol(:) = 1.e6
996        max_resol(:) = -1.
997        !---
998        !--- All grid variables
999        !---
1000        !-
1001        !- Get variables to help describe the grid
1002        CALL flinquery_var(force_id, 'neighboursNN', neighbours_exists)
1003        IF ( .NOT. neighbours_exists ) THEN
1004           CALL ipslerr_p (3,'forcing_read_interpol', &
1005 &          'Could get neighbours in a watchout file :',TRIM(filename), &
1006 &          '(Problem with file ?)')
1007        ELSE
1008           WRITE(numout,*) "Watchout file contains neighbours and resolutions."
1009        ENDIF
1010        !
1011        fneighbours(:,:,:) = undef_int
1012        !
1013        !- once the points are gathered.
1014        CALL flinget_buffer (force_id,'neighboursNN',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1015        CALL forcing_zoom(data_full, tmpdata)
1016        WHERE(tmpdata(:,:) < undef_int)
1017           fneighbours(:,:,1) = NINT(tmpdata(:,:))
1018        ENDWHERE
1019        !
1020        CALL flinget_buffer (force_id,'neighboursNE',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1021        CALL forcing_zoom(data_full, tmpdata)
1022        WHERE(tmpdata(:,:) < undef_int)
1023           fneighbours(:,:,2) = NINT(tmpdata(:,:))
1024        ENDWHERE
1025        !
1026        CALL flinget_buffer (force_id,'neighboursEE',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1027        CALL forcing_zoom(data_full, tmpdata)
1028        WHERE(tmpdata(:,:) < undef_int)
1029           fneighbours(:,:,3) = NINT(tmpdata(:,:))
1030        ENDWHERE
1031        !
1032        CALL flinget_buffer (force_id,'neighboursSE',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1033        CALL forcing_zoom(data_full, tmpdata)
1034        WHERE(tmpdata(:,:) < undef_int)
1035           fneighbours(:,:,4) = NINT(tmpdata(:,:))
1036        ENDWHERE
1037        !
1038        CALL flinget_buffer (force_id,'neighboursSS',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1039        CALL forcing_zoom(data_full, tmpdata)
1040        WHERE(tmpdata(:,:) < undef_int)
1041           fneighbours(:,:,5) = NINT(tmpdata(:,:))
1042        ENDWHERE
1043        !
1044        CALL flinget_buffer (force_id,'neighboursSW',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1045        CALL forcing_zoom(data_full, tmpdata)
1046        WHERE(tmpdata(:,:) < undef_int)
1047           fneighbours(:,:,6) = NINT(tmpdata(:,:))
1048        ENDWHERE
1049        !
1050        CALL flinget_buffer (force_id,'neighboursWW',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1051        CALL forcing_zoom(data_full, tmpdata)
1052        WHERE(tmpdata(:,:) < undef_int)
1053           fneighbours(:,:,7) = NINT(tmpdata(:,:))
1054        ENDWHERE
1055        !
1056        CALL flinget_buffer (force_id,'neighboursNW',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1057        CALL forcing_zoom(data_full, tmpdata)
1058        WHERE(tmpdata(:,:) < undef_int)
1059           fneighbours(:,:,8) = NINT(tmpdata(:,:))
1060        ENDWHERE
1061        !
1062        ! now, resolution of the grid
1063        CALL flinget_buffer (force_id,'resolutionX',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1064        CALL forcing_zoom(data_full, tmpdata)
1065        fresolution(:,:,1) = tmpdata(:,:)
1066        !
1067        CALL flinget_buffer (force_id,'resolutionY',iim_full, jjm_full, llm_full, ttm,  1, 1, data_full)
1068        CALL forcing_zoom(data_full, tmpdata)
1069        fresolution(:,:,2) = tmpdata(:,:)
1070        !
1071        DO ik=1,nbindex
1072           !
1073           j = ((kindex(ik)-1)/iim) + 1
1074           i = (kindex(ik) - (j-1)*iim)
1075           !-
1076           !- Store variable to help describe the grid
1077           !- once the points are gathered.
1078           !-
1079           neighbours(ik,:) = fneighbours(i,j,:) 
1080           !
1081           resolution(ik,:) = fresolution(i,j,:)
1082           !
1083       
1084        ENDDO
1085        CALL gather(neighbours,neighbours_g)
1086        CALL gather(resolution,resolution_g)
1087        min_resol(1) = MINVAL(resolution(:,1))
1088        min_resol(2) = MAXVAL(resolution(:,2))
1089        max_resol(1) = MAXVAL(resolution(:,1))
1090        max_resol(2) = MAXVAL(resolution(:,2))
1091        !
1092        area(:) = resolution(:,1)*resolution(:,2)
1093        CALL gather(area,area_g)
1094!--
1095        DEALLOCATE (tmpdata)
1096     ENDIF
1097     WRITE(numout,*) 'contfrac : ', MINVAL(fcontfrac), MAXVAL(fcontfrac)
1098!---
1099   ENDIF
1100!---
1101   IF (check) THEN
1102      WRITE(numout,*) &
1103           & 'The dates : ',itau_read,itau_split,itau_read_nm1,itau_read_n
1104   ENDIF
1105!---
1106!--- Here we do the work in case only interpolation is needed.
1107!---
1108   IF ( initialized .AND. interpol ) THEN
1109!---
1110      IF ( daily_interpol ) THEN
1111
1112         IF (split > 1) THEN
1113            IF ( itau_split <= (split/2.) ) THEN
1114               rw = REAL(itau_split+split/2.)/split
1115            ELSE
1116               rw = REAL(itau_split-split/2.)/split
1117            ENDIF
1118         ELSE
1119            rw = 1.
1120         ENDIF
1121
1122         IF ((last_read == 0) .OR. ( rw==(1./split)) ) THEN
1123   !---
1124   !-----   Start or Restart
1125            IF (last_read == 0) THEN
1126               ! Case of a restart or a shift in the forcing file.
1127               IF (itau_read > 1) THEN
1128                  itau_read_nm1=itau_read-1
1129                  CALL forcing_just_read (iim, jjm, zlev_nm1, zlevuv_nm1, ttm, itau_read_nm1, itau_read_nm1, &
1130                       &  swdown_nm1, rainf_nm1, snowf_nm1, tmin_nm1, &
1131                       &  u_nm1, v_nm1, qair_nm1, pb_nm1, lwdown_nm1, &
1132                       &  SWnet_nm1, Eair_nm1, petAcoef_nm1, peqAcoef_nm1, petBcoef_nm1, peqBcoef_nm1, cdrag_nm1, ccanopy_nm1, &
1133                       &  force_id, wind_N_exists, check)
1134                  CALL forcing_just_read_tmax (iim, jjm, ttm, itau_read_nm1, itau_read_nm1, tmax_nm1, force_id )
1135               ! Case of a simple start.
1136               ELSE
1137                  itau_read_nm1 = un
1138                  WRITE(numout,*) "we will use the forcing of the first day to initialize "
1139                  CALL forcing_just_read (iim, jjm, zlev_nm1, zlevuv_nm1, ttm, itau_read_nm1, itau_read_nm1, &
1140                       &  swdown_nm1, rainf_nm1, snowf_nm1, tmin_nm1, &
1141                       &  u_nm1, v_nm1, qair_nm1, pb_nm1, lwdown_nm1, &
1142                       &  SWnet_nm1, Eair_nm1, petAcoef_nm1, peqAcoef_nm1, petBcoef_nm1, peqBcoef_nm1, cdrag_nm1, ccanopy_nm1, &
1143                       &  force_id, wind_N_exists, check)
1144                  CALL forcing_just_read_tmax (iim, jjm, ttm, itau_read_nm1, itau_read_nm1, tmax_nm1, force_id )
1145               ENDIF
1146               tmin_nm2(:,:)=tmin_nm1(:,:)
1147               IF ( dt_force .GT. 3600. ) THEN
1148                  mean_coszang(:,:) = 0.0
1149                  daylength_n(:,:) = 0.
1150                  DO is=1,split
1151                     !MM we compute mean SWdown between t and t+Dt then I take t+Dt/2.
1152                     julian = julian_for+((is-0.5)/split)*dt_force/one_day
1153      !!$               julian = julian_for+(FLOAT(is)/split)*dt_force/one_day
1154                     CALL solarang (julian, julian0, iim, jjm, lon*0, lat, coszang)
1155                     mean_coszang(:,:) = mean_coszang(:,:)+coszang(:,:)
1156                     WHERE( coszang(:,:) > 0. ) 
1157                        daylength_n(:,:)=daylength_n(:,:)+1./split*24
1158                     ENDWHERE
1159                  ENDDO
1160                  mean_coszang(:,:) = mean_coszang(:,:)/split
1161                  daylength_nm1(:,:)=daylength_n(:,:)
1162      !            WRITE(*,*) "mean_coszang =",MAXVAL(mean_coszang)
1163               ENDIF
1164            ELSE
1165   !-----   Normal mode : copy old step
1166               swdown_nm1(:,:) = swdown_n(:,:)
1167               rainf_nm1(:,:) = rainf_n(:,:)
1168               snowf_nm1(:,:)  = snowf_n(:,:) 
1169               tair_nm1(:,:)   = tair_n(:,:)
1170               u_nm1(:,:)      = u_n(:,:)
1171               v_nm1(:,:)      = v_n(:,:)
1172               qair_nm1(:,:)   = qair_n(:,:)
1173               pb_nm1(:,:)     = pb_n(:,:)
1174               lwdown_nm1(:,:) = lwdown_n(:,:)
1175               tmin_nm2(:,:)   = tmin_nm1(:,:)
1176               tmin_nm1(:,:)   = tmin_n(:,:)
1177               tmax_nm1(:,:)   = tmax_n(:,:)
1178
1179               IF (is_watchout) THEN
1180                  zlev_nm1(:,:)   = zlev_n(:,:)
1181                  zlevuv_nm1(:,:) = zlevuv_n(:,:)
1182                  ! Net surface short-wave flux
1183                  SWnet_nm1(:,:) = SWnet_n(:,:)
1184                  ! Air potential energy
1185                  Eair_nm1(:,:)   = Eair_n(:,:)
1186                  ! Coeficients A from the PBL resolution for T
1187                  petAcoef_nm1(:,:) = petAcoef_n(:,:)
1188                  ! Coeficients A from the PBL resolution for q
1189                  peqAcoef_nm1(:,:) = peqAcoef_n(:,:)
1190                  ! Coeficients B from the PBL resolution for T
1191                  petBcoef_nm1(:,:) = petBcoef_n(:,:)
1192                  ! Coeficients B from the PBL resolution for q
1193                  peqBcoef_nm1(:,:) = peqBcoef_n(:,:)
1194                  ! Surface drag
1195                  cdrag_nm1(:,:) = cdrag_n(:,:)
1196                  ! CO2 concentration in the canopy
1197                  ccanopy_nm1(:,:) = ccanopy_n(:,:)
1198               ENDIF
1199               itau_read_nm1 = itau_read_n
1200            ENDIF
1201   !-----
1202   !-----
1203            IF(last_read==0)THEN
1204               itau_read_n = itau_read
1205            ELSE
1206               itau_read_n = itau_read+1
1207            ENDIF
1208
1209            IF (itau_read_n > ttm) THEN
1210               WRITE(numout,*) 'WARNING --WARNING --WARNING --WARNING '
1211               WRITE(numout,*) &
1212                    &  'WARNING : We are going back to the start of the file'
1213               itau_read_n =1
1214            ENDIF
1215            IF (check) THEN
1216               WRITE(numout,*) &
1217                    & 'The dates 2 : ',itau_read,itau_split,itau_read_nm1,itau_read_n
1218            ENDIF
1219   !-----
1220   !----- Get a reduced julian day !
1221   !----- This is needed because we lack the precision on 32 bit machines.
1222   !-----
1223            IF ( dt_force .GT. 3600. ) THEN
1224               julian_for = itau2date(itau_read-1, date0, dt_force)
1225               CALL ju2ymds (julian_for, yy, mm, dd, ss)
1226   
1227               ! first day of this year
1228               CALL ymds2ju (yy,1,1,0.0, julian0)
1229   !-----
1230               IF (check) THEN
1231                  WRITE(numout,*) 'Forcing for Julian day ',julian_for,'is read'
1232                  WRITE(numout,*) 'Date for this day ',yy,' / ',mm,' / ',dd,"  ",ss
1233               ENDIF
1234            ENDIF
1235   !-----
1236            CALL forcing_just_read (iim, jjm, zlev_n, zlevuv_n, ttm, itau_read_n, itau_read_n, &
1237                 &  swdown_n, rainf_n, snowf_n, tmin_n, &
1238                 &  u_n, v_n, qair_n, pb_n, lwdown_n, &
1239                 &  SWnet_n, Eair_n, petAcoef_n, peqAcoef_n, petBcoef_n, peqBcoef_n, cdrag_n, ccanopy_n, &
1240                 &  force_id, wind_N_exists, check)
1241            CALL forcing_just_read_tmax (iim, jjm, ttm, itau_read_n, itau_read_n, tmax_n, force_id )
1242
1243   !---
1244            last_read = itau_read_n
1245   !-----
1246   !----- Compute mean solar angle for the comming period
1247   !-----
1248            IF (check) WRITE(numout,*) 'Going into  solarang', split, one_day
1249   !-----
1250
1251   !-----
1252         ENDIF
1253   !---
1254         IF ( itau_split == 1. ) THEN
1255            IF ( dt_force .GT. 3600. ) THEN
1256               mean_coszang(:,:) = 0.0
1257               daylength_nm1(:,:)=daylength_n(:,:)
1258               daylength_n(:,:) = 0.
1259               DO is=1,split
1260                  !MM we compute mean SWdown between t and t+Dt then I take t+Dt/2.
1261                  julian = julian_for+((is-0.5)/split)*dt_force/one_day
1262   !!$               julian = julian_for+(FLOAT(is)/split)*dt_force/one_day
1263                  CALL solarang (julian, julian0, iim, jjm, lon*0, lat, coszang)
1264                  mean_coszang(:,:) = mean_coszang(:,:)+coszang(:,:)
1265                  WHERE( coszang(:,:) > 0. ) 
1266                     daylength_n(:,:)=daylength_n(:,:)+1./split*24
1267                  ENDWHERE
1268               ENDDO
1269               mean_coszang(:,:) = mean_coszang(:,:)/split
1270   !            WRITE(*,*) "mean_coszang =",MAXVAL(mean_coszang)
1271            ENDIF
1272         ENDIF
1273 
1274   !--- Do the interpolation
1275         IF (check) WRITE(numout,*) 'Doing the interpolation between time steps'
1276   !---
1277
1278         IF (check) WRITE(numout,*) 'Coeff of interpollation : ',rw
1279   !---
1280
1281         pb(:,:) = (pb_n(:,:)-pb_nm1(:,:))*rw + pb_nm1(:,:)
1282         u(:,:)  = (u_n(:,:)-u_nm1(:,:))*rw + u_nm1(:,:)
1283         v(:,:)  = (v_n(:,:)-v_nm1(:,:))*rw + v_nm1(:,:)
1284
1285   !--- Take care of the height of the vertical levels
1286         zlev(:,:) = (zlev_n(:,:)-zlev_nm1(:,:))*rw + zlev_nm1(:,:) 
1287         zlevuv(:,:) = (zlevuv_n(:,:)-zlevuv_nm1(:,:))*rw + zlevuv_nm1(:,:) 
1288
1289         hour=REAL(itau_split)/split*24
1290         startday_n(:,:)=12.-daylength_n(:,:)/2.
1291         startday_nm1(:,:)=12.-daylength_nm1(:,:)/2.
1292
1293         WHERE ( ( hour >= startday_n(:,:) ) .AND. ( hour > 12) .AND. ( hour <= 14) )
1294            tair(:,:)=(tmax_nm1(:,:)-tmin_nm1(:,:))/2 * ( sin(pi/(14-startday_n(:,:))*(hour-0.5* &
1295           &    (14.-startday_n(:,:))-startday_n(:,:))) )+ (tmax_nm1(:,:)+tmin_nm1(:,:))/2.
1296         ELSEWHERE( ( hour >= startday_n(:,:) ) .AND. ( hour <= 12) )
1297            tair(:,:)=(tmax_n(:,:)-tmin_n(:,:))/2 * ( sin(pi/(14-startday_n(:,:))*(hour-0.5* &
1298           &    (14.-startday_n(:,:))-startday_n(:,:))) )+ (tmax_n(:,:)+tmin_n(:,:))/2.
1299         ELSEWHERE ( hour < startday_n(:,:) )
1300            tair(:,:)=(tmax_nm1(:,:)-tmin_n(:,:))/2.*sin(pi/(24.-14.+startday_nm1(:,:) )* &
1301           &    (hour + 24.+0.5*(24.-14.+startday_nm1(:,:) )-14.))+(tmax_nm1(:,:)+tmin_n(:,:))/2.
1302         ELSEWHERE
1303            tair(:,:)=(tmax_nm1(:,:)-tmin_n(:,:))/2.*sin(pi/(24.-14.+startday_n(:,:))*(hour+0.5* &
1304           &    (24.-14.+startday_n(:,:))-14.))+(tmax_nm1(:,:)+tmin_n(:,:))/2.
1305         ENDWHERE
1306
1307         CALL weathgen_qsat_2d (iim,jjm,tmin_n,pb,qsattmin_n)
1308         CALL weathgen_qsat_2d (iim,jjm,tmin_nm1,pb,qsattmin_nm1)
1309         CALL weathgen_qsat_2d (iim,jjm,tair,pb,qsatta)
1310
1311         !---
1312         qmin_nm1(:,:) = MIN(qair_nm1(:,:),0.99*qsattmin_nm1(:,:))
1313         qmin_n(:,:) = MIN(qair_n(:,:),0.99*qsattmin_n(:,:))
1314         qmax_nm1(:,:) = (qair_nm1(:,:)-qmin_nm1(:,:)) + qair_nm1(:,:)
1315         qmax_n(:,:) = (qair_n(:,:)-qmin_n(:,:)) + qair_n(:,:)
1316
1317         qsa(:,:)  = 0.99*qsatta(:,:)
1318
1319
1320         WHERE ( ( hour >= startday_n(:,:) ) .AND. ( hour > 12) .AND. ( hour <= 14) )
1321            qair(:,:)=MIN(qsa(:,:),(qmax_nm1(:,:)-qmin_nm1(:,:))/2 * ( sin(pi/(14-startday_n(:,:))*(hour-0.5* &
1322           &    (14.-startday_n(:,:))-startday_n(:,:))) )+ (qmax_nm1(:,:)+qmin_nm1(:,:))/2.)
1323         ELSEWHERE( ( hour >= startday_n(:,:) ) .AND. ( hour <= 12) )
1324            qair(:,:)=MIN(qsa(:,:),(qmax_n(:,:)-qmin_n(:,:))/2 * ( sin(pi/(14-startday_n(:,:))*(hour-0.5* &
1325           &    (14.-startday_n(:,:))-startday_n(:,:))) )+ (qmax_n(:,:)+qmin_n(:,:))/2.)
1326         ELSEWHERE ( hour < startday_n(:,:) )
1327            qair(:,:)=MIN(qsa(:,:),(qmax_nm1(:,:)-qmin_n(:,:))/2.*sin(pi/(24.-14.+startday_nm1(:,:) )* &
1328           &    (hour + 24.+0.5*(24.-14.+startday_nm1(:,:) )-14.))+(qmax_nm1(:,:)+qmin_n(:,:))/2.)
1329         ELSEWHERE
1330            qair(:,:)=MIN(qsa(:,:),(qmax_nm1(:,:)-qmin_n(:,:))/2.*sin(pi/(24.-14.+startday_n(:,:))*(hour+0.5* &
1331           &    (24.-14.+startday_n(:,:))-14.))+(qmax_nm1(:,:)+qmin_n(:,:))/2.)
1332         ENDWHERE
1333
1334         IF (is_watchout) THEN
1335            SWnet(:,:) = (SWnet_n(:,:)-SWnet_nm1(:,:))*rw + SWnet_nm1(:,:)
1336            Eair(:,:) = (Eair_n(:,:)-Eair_nm1(:,:))*rw + Eair_nm1(:,:)
1337            petAcoef(:,:) = (petAcoef_n(:,:)-petAcoef_nm1(:,:))*rw + petAcoef_nm1(:,:)
1338            peqAcoef(:,:) = (peqAcoef_n(:,:)-peqAcoef_nm1(:,:))*rw + peqAcoef_nm1(:,:)
1339            petBcoef(:,:) = (petBcoef_n(:,:)-petBcoef_nm1(:,:))*rw + petBcoef_nm1(:,:)
1340            peqBcoef(:,:) = (peqBcoef_n(:,:)-peqBcoef_nm1(:,:))*rw + peqBcoef_nm1(:,:)
1341            cdrag(:,:) = (cdrag_n(:,:)-cdrag_nm1(:,:))*rw + cdrag_nm1(:,:)
1342            ccanopy(:,:) = (ccanopy_n(:,:)-ccanopy_nm1(:,:))*rw + ccanopy_nm1(:,:)
1343         ENDIF
1344   !---
1345   !--- Here we need to allow for an option
1346   !--- where radiative energy is conserved
1347   !---
1348         IF ( netrad_cons ) THEN
1349            lwdown(:,:) = lwdown_n(:,:)
1350         ELSE
1351            lwdown(:,:) = (lwdown_n(:,:)-lwdown_nm1(:,:))*rw + lwdown_nm1(:,:)
1352         ENDIF
1353   !---
1354   !--- For the solar radiation we decompose the mean value
1355   !--- using the zenith angle of the sun if the time step in the forcing data is
1356   !---- more than an hour. Else we use the standard linera interpolation
1357   !----
1358         IF (check) WRITE(numout,*) 'Ready to deal with the solar radiation'
1359   !----
1360         IF ( dt_force .GT. 3600. ) THEN
1361   !---
1362            IF ( netrad_cons ) THEN
1363               WRITE(numout,*) 'Solar radiation can not be conserved with a timestep of ', dt_force
1364            ENDIF
1365   !---
1366            !MM we compute mean SWdown between t and t+Dt then I take t+Dt/2.
1367            julian = julian_for + (itau_split-0.5)/split*dt_force/one_day
1368   !!$         julian = julian_for + rw*dt_force/one_day
1369            IF (check) THEN
1370               WRITE(numout,'(a,f20.10,2I3)') &
1371                    &  'JULIAN BEFORE SOLARANG : ',julian,itau_split,split
1372            ENDIF
1373   !---
1374            CALL solarang(julian, julian0, iim, jjm, lon*0, lat, coszang)
1375   !---
1376           
1377            WHERE ((mean_coszang(:,:) > 0.) .AND. (hour <= 12 ))
1378               swdown(:,:) = swdown_n(:,:) *coszang(:,:)/mean_coszang(:,:)
1379            ELSEWHERE ((mean_coszang(:,:) > 0.) .AND. (hour > 12 ))
1380               swdown(:,:) = swdown_nm1(:,:) *coszang(:,:)/mean_coszang(:,:)
1381            ELSEWHERE
1382               swdown(:,:) = 0.0
1383            END WHERE
1384   !---
1385            WHERE (swdown(:,:) > 2000. )
1386               swdown(:,:) = 2000.
1387            END WHERE
1388   !---
1389         ELSE
1390   !---
1391            IF ( netrad_cons ) THEN
1392               swdown(:,:) = swdown_n(:,:)
1393            ELSE
1394               swdown(:,:) = (swdown_n(:,:)-swdown_nm1(:,:))*rw + swdown_nm1(:,:)
1395            ENDIF
1396   !---
1397         ENDIF
1398   !---
1399         IF (check) THEN
1400            WRITE(numout,*) '__ Forcing read at ',itau_split,' :',i_test, j_test
1401            WRITE(numout,*) 'SWdown  : ',swdown_nm1(i_test, j_test), &
1402                 &           ' < ', swdown(i_test, j_test), ' < ', swdown_n(i_test, j_test)
1403            IF (is_watchout) THEN
1404               WRITE(numout,*) 'SWnet  : ',swnet_nm1(i_test, j_test), &
1405                    &           ' < ', swnet(i_test, j_test), ' < ', swnet_n(i_test, j_test)
1406               WRITE(numout,*) 'levels  :',zlev_nm1(i_test, j_test), &
1407                    &           ' < ', zlev(i_test, j_test), ' < ', zlev_n(i_test, j_test)
1408               WRITE(numout,*) 'EAIR  :',Eair_nm1(i_test, j_test), &
1409                    &           ' < ', eair(i_test, j_test), ' < ', Eair_n(i_test, j_test)
1410            ENDIF
1411            WRITE(numout,*) 'TAIR  :',tair_nm1(i_test, j_test), &
1412                 &           ' < ', tair(i_test, j_test), ' < ', tair_n(i_test, j_test)
1413            WRITE(numout,*) 'QAIR  :',qair_nm1(i_test, j_test), &
1414                 &           ' < ', qair(i_test, j_test), ' < ', qair_n(i_test, j_test)
1415            WRITE(numout,*) 'U  :',u_nm1(i_test, j_test), &
1416                 &           ' < ', u(i_test, j_test), ' < ', u_n(i_test, j_test)
1417            WRITE(numout,*) 'V  :',v_nm1(i_test, j_test), &
1418                 &           ' < ', v(i_test, j_test), ' < ', v_n(i_test, j_test)
1419         ENDIF
1420   !---
1421   !--- For precip we suppose that the rain
1422   !--- is the sum over the next 6 hours
1423   !---
1424         WHERE ((itau_split <= nb_spread).AND.(hour<=12).AND.(tair(:,:)>=273.15)) 
1425            rainf(:,:) = rainf_n(:,:) *(split/REAL(nb_spread))
1426            snowf(:,:) = 0.0
1427         ELSEWHERE ((itau_split <= nb_spread).AND.(hour<=12).AND.(tair(:,:)<273.15)) 
1428            snowf(:,:) = rainf_n(:,:) *(split/REAL(nb_spread))
1429            rainf(:,:) = 0.0
1430         ELSEWHERE ((itau_split <= nb_spread).AND.(hour>12).AND.(tair(:,:)>=273.15)) 
1431            rainf(:,:) = rainf_nm1(:,:) *(split/REAL(nb_spread))
1432            snowf(:,:) = 0.0
1433         ELSEWHERE ((itau_split <= nb_spread).AND.(hour>12).AND.(tair(:,:)<273.15)) 
1434            snowf(:,:) = rainf_nm1(:,:) *(split/REAL(nb_spread))
1435            rainf(:,:) = 0.0
1436         ELSEWHERE
1437            snowf(:,:) = 0.0
1438            rainf(:,:) = 0.0
1439         ENDWHERE
1440
1441         IF (check) THEN
1442            WRITE(numout,*) '__ Forcing read at ',itau_split,' :'
1443            WRITE(numout,*) 'Rainf  : ',rainf_nm1(i_test, j_test), &
1444                 &           ' < ', rainf(i_test, j_test), ' < ', rainf_n(i_test, j_test)
1445            WRITE(numout,*) 'Snowf  : ',snowf_nm1(i_test, j_test), &
1446                 &           ' < ', snowf(i_test, j_test), ' < ', snowf_n(i_test, j_test)
1447         ENDIF
1448   !---
1449
1450
1451      ELSE
1452     
1453         IF (itau_read /= last_read) THEN
1454   !---
1455   !-----   Start or Restart
1456            IF (itau_read_n == 0) THEN
1457               ! Case of a restart or a shift in the forcing file.
1458               IF (itau_read > 1) THEN
1459                  itau_read_nm1=itau_read-1
1460                  CALL forcing_just_read (iim, jjm, zlev_nm1, zlevuv_nm1, ttm, itau_read_nm1, itau_read_nm1, &
1461                       &  swdown_nm1, rainf_nm1, snowf_nm1, tair_nm1, &
1462                       &  u_nm1, v_nm1, qair_nm1, pb_nm1, lwdown_nm1, &
1463                       &  SWnet_nm1, Eair_nm1, petAcoef_nm1, peqAcoef_nm1, petBcoef_nm1, peqBcoef_nm1, cdrag_nm1, ccanopy_nm1, &
1464                       &  force_id, wind_N_exists, check)
1465               ! Case of a simple start.
1466               ELSE IF (dt_force*ttm > one_day-1. ) THEN
1467                  ! if the forcing file contains at least 24 hours,
1468                  ! we will use the last forcing step of the first day
1469                  ! as initiale condition to prevent first shift off reading.
1470                  itau_read_nm1 = NINT (one_day/dt_force)
1471                  WRITE(numout,*) "the forcing file contains 24 hours :",dt_force*ttm,one_day-1.
1472                  WRITE(numout,*) "we will use the last forcing step of the first day : itau_read_nm1 ",itau_read_nm1
1473                  CALL forcing_just_read (iim, jjm, zlev_nm1, zlevuv_nm1, ttm, itau_read_nm1, itau_read_nm1, &
1474                       &  swdown_nm1, rainf_nm1, snowf_nm1, tair_nm1, &
1475                       &  u_nm1, v_nm1, qair_nm1, pb_nm1, lwdown_nm1, &
1476                       &  SWnet_nm1, Eair_nm1, petAcoef_nm1, peqAcoef_nm1, petBcoef_nm1, peqBcoef_nm1, cdrag_nm1, ccanopy_nm1, &
1477                       &  force_id, wind_N_exists, check)
1478               ELSE
1479                  ! if the forcing file contains less than 24 hours,
1480                  ! just say error !
1481                  CALL ipslerr_p(3,'forcing_read_interpol', &
1482      &             'The forcing file contains less than 24 hours !', &
1483      &             'We can''t intialize interpolation with such a file.','') 
1484               ENDIF
1485            ELSE
1486   !-----   Normal mode : copy old step
1487               swdown_nm1(:,:) = swdown_n(:,:)
1488               rainf_nm1(:,:) = rainf_n(:,:)
1489               snowf_nm1(:,:)  = snowf_n(:,:) 
1490               tair_nm1(:,:)   = tair_n(:,:)
1491               u_nm1(:,:)      = u_n(:,:)
1492               v_nm1(:,:)      = v_n(:,:)
1493               qair_nm1(:,:)   = qair_n(:,:)
1494               pb_nm1(:,:)     = pb_n(:,:)
1495               lwdown_nm1(:,:) = lwdown_n(:,:)
1496               IF (is_watchout) THEN
1497                  zlev_nm1(:,:)   = zlev_n(:,:)
1498                  ! Net surface short-wave flux
1499                  SWnet_nm1(:,:) = SWnet_n(:,:)
1500                  ! Air potential energy
1501                  Eair_nm1(:,:)   = Eair_n(:,:)
1502                  ! Coeficients A from the PBL resolution for T
1503                  petAcoef_nm1(:,:) = petAcoef_n(:,:)
1504                  ! Coeficients A from the PBL resolution for q
1505                  peqAcoef_nm1(:,:) = peqAcoef_n(:,:)
1506                  ! Coeficients B from the PBL resolution for T
1507                  petBcoef_nm1(:,:) = petBcoef_n(:,:)
1508                  ! Coeficients B from the PBL resolution for q
1509                  peqBcoef_nm1(:,:) = peqBcoef_n(:,:)
1510                  ! Surface drag
1511                  cdrag_nm1(:,:) = cdrag_n(:,:)
1512                  ! CO2 concentration in the canopy
1513                  ccanopy_nm1(:,:) = ccanopy_n(:,:)
1514               ENDIF
1515               itau_read_nm1 = itau_read_n
1516            ENDIF
1517   !-----
1518            itau_read_n = itau_read
1519            IF (itau_read_n > ttm) THEN
1520               WRITE(numout,*) 'WARNING --WARNING --WARNING --WARNING '
1521               WRITE(numout,*) &
1522                    &  'WARNING : We are going back to the start of the file'
1523               itau_read_n =1
1524            ENDIF
1525            IF (check) THEN
1526               WRITE(numout,*) &
1527                    & 'The dates 2 : ',itau_read,itau_split,itau_read_nm1,itau_read_n
1528            ENDIF
1529   !-----
1530   !----- Get a reduced julian day !
1531   !----- This is needed because we lack the precision on 32 bit machines.
1532   !-----
1533            IF ( dt_force .GT. 3600. ) THEN
1534               julian_for = itau2date(itau_read-1, date0, dt_force)
1535               CALL ju2ymds (julian_for, yy, mm, dd, ss)
1536   
1537               ! first day of this year
1538               CALL ymds2ju (yy,1,1,0.0, julian0)
1539   !-----
1540               IF (check) THEN
1541                  WRITE(numout,*) 'Forcing for Julian day ',julian_for,'is read'
1542                  WRITE(numout,*) 'Date for this day ',yy,' / ',mm,' / ',dd,"  ",ss
1543               ENDIF
1544            ENDIF
1545   !-----
1546            CALL forcing_just_read (iim, jjm, zlev_n, zlevuv_n, ttm, itau_read_n, itau_read_n, &
1547                 &  swdown_n, rainf_n, snowf_n, tair_n, &
1548                 &  u_n, v_n, qair_n, pb_n, lwdown_n, &
1549                 &  SWnet_n, Eair_n, petAcoef_n, peqAcoef_n, petBcoef_n, peqBcoef_n, cdrag_n, ccanopy_n, &
1550                 &  force_id, wind_N_exists, check)
1551   !---
1552            last_read = itau_read_n
1553   !-----
1554   !----- Compute mean solar angle for the comming period
1555   !-----
1556            IF (check) WRITE(numout,*) 'Going into  solarang', split, one_day
1557   !-----
1558            IF ( dt_force .GT. 3600. ) THEN
1559               mean_coszang(:,:) = 0.0
1560               DO is=1,split
1561                  !MM we compute mean SWdown between t and t+Dt then I take t+Dt/2.
1562                  julian = julian_for+((is-0.5)/split)*dt_force/one_day
1563   !!$               julian = julian_for+(FLOAT(is)/split)*dt_force/one_day
1564                  CALL solarang (julian, julian0, iim, jjm, lon, lat, coszang)
1565                  mean_coszang(:,:) = mean_coszang(:,:)+coszang(:,:)
1566               ENDDO
1567               mean_coszang(:,:) = mean_coszang(:,:)/split
1568   !            WRITE(*,*) "mean_coszang =",MAXVAL(mean_coszang)
1569            ENDIF
1570   !-----
1571         ENDIF
1572   !---
1573   !--- Do the interpolation
1574         IF (check) WRITE(numout,*) 'Doing the interpolation between time steps'
1575   !---
1576         IF (split > 1) THEN
1577            rw = REAL(itau_split)/split
1578         ELSE
1579            rw = 1.
1580         ENDIF
1581         IF (check) WRITE(numout,*) 'Coeff of interpollation : ',rw
1582   !---
1583         qair(:,:) = (qair_n(:,:)-qair_nm1(:,:))*rw + qair_nm1(:,:)
1584         tair(:,:) = (tair_n(:,:)-tair_nm1(:,:))*rw + tair_nm1(:,:)
1585         pb(:,:) = (pb_n(:,:)-pb_nm1(:,:))*rw + pb_nm1(:,:)
1586         u(:,:)  = (u_n(:,:)-u_nm1(:,:))*rw + u_nm1(:,:)
1587         v(:,:)  = (v_n(:,:)-v_nm1(:,:))*rw + v_nm1(:,:)
1588         IF (is_watchout) THEN
1589            zlev(:,:) = (zlev_n(:,:)-zlev_nm1(:,:))*rw + zlev_nm1(:,:)
1590            zlevuv(:,:) = zlev(:,:)
1591            SWnet(:,:) = (SWnet_n(:,:)-SWnet_nm1(:,:))*rw + SWnet_nm1(:,:)
1592            Eair(:,:) = (Eair_n(:,:)-Eair_nm1(:,:))*rw + Eair_nm1(:,:)
1593            petAcoef(:,:) = (petAcoef_n(:,:)-petAcoef_nm1(:,:))*rw + petAcoef_nm1(:,:)
1594            peqAcoef(:,:) = (peqAcoef_n(:,:)-peqAcoef_nm1(:,:))*rw + peqAcoef_nm1(:,:)
1595            petBcoef(:,:) = (petBcoef_n(:,:)-petBcoef_nm1(:,:))*rw + petBcoef_nm1(:,:)
1596            peqBcoef(:,:) = (peqBcoef_n(:,:)-peqBcoef_nm1(:,:))*rw + peqBcoef_nm1(:,:)
1597            cdrag(:,:) = (cdrag_n(:,:)-cdrag_nm1(:,:))*rw + cdrag_nm1(:,:)
1598            ccanopy(:,:) = (ccanopy_n(:,:)-ccanopy_nm1(:,:))*rw + ccanopy_nm1(:,:)
1599         ENDIF
1600   !---
1601   !--- Here we need to allow for an option
1602   !--- where radiative energy is conserved
1603   !---
1604         IF ( netrad_cons ) THEN
1605            lwdown(:,:) = lwdown_n(:,:)
1606         ELSE
1607            lwdown(:,:) = (lwdown_n(:,:)-lwdown_nm1(:,:))*rw + lwdown_nm1(:,:)
1608         ENDIF
1609   !---
1610   !--- For the solar radiation we decompose the mean value
1611   !--- using the zenith angle of the sun if the time step in the forcing data is
1612   !---- more than an hour. Else we use the standard linera interpolation
1613   !----
1614         IF (check) WRITE(numout,*) 'Ready to deal with the solar radiation'
1615   !----
1616         IF ( dt_force .GT. 3600. ) THEN
1617   !---
1618            IF ( netrad_cons ) THEN
1619               WRITE(numout,*) 'Solar radiation can not be conserved with a timestep of ', dt_force
1620            ENDIF
1621   !---
1622            !MM we compute mean SWdown between t and t+Dt then I take t+Dt/2.
1623            julian = julian_for + (itau_split-0.5)/split*dt_force/one_day
1624   !!$         julian = julian_for + rw*dt_force/one_day
1625            IF (check) THEN
1626               WRITE(numout,'(a,f20.10,2I3)') &
1627                    &  'JULIAN BEFORE SOLARANG : ',julian,itau_split,split
1628            ENDIF
1629   !---
1630            CALL solarang(julian, julian0, iim, jjm, lon, lat, coszang)
1631   !---
1632            WHERE (mean_coszang(:,:) > 0.)
1633               swdown(:,:) = swdown_n(:,:) *coszang(:,:)/mean_coszang(:,:)
1634            ELSEWHERE
1635               swdown(:,:) = 0.0
1636            END WHERE
1637   !---
1638            WHERE (swdown(:,:) > 2000. )
1639               swdown(:,:) = 2000.
1640            END WHERE
1641   !---
1642         ELSE
1643   !---
1644            IF ( netrad_cons ) THEN
1645               swdown(:,:) = swdown_n(:,:)
1646            ELSE
1647               swdown(:,:) = (swdown_n(:,:)-swdown_nm1(:,:))*rw + swdown_nm1(:,:)
1648            ENDIF
1649   !---
1650         ENDIF
1651   !---
1652         IF (check) THEN
1653            WRITE(numout,*) '__ Forcing read at ',itau_split,' :',i_test, j_test
1654            WRITE(numout,*) 'SWdown  : ',swdown_nm1(i_test, j_test), &
1655                 &           ' < ', swdown(i_test, j_test), ' < ', swdown_n(i_test, j_test)
1656            IF (is_watchout) THEN
1657               WRITE(numout,*) 'SWnet  : ',swnet_nm1(i_test, j_test), &
1658                    &           ' < ', swnet(i_test, j_test), ' < ', swnet_n(i_test, j_test)
1659               WRITE(numout,*) 'levels  :',zlev_nm1(i_test, j_test), &
1660                    &           ' < ', zlev(i_test, j_test), ' < ', zlev_n(i_test, j_test)
1661               WRITE(numout,*) 'EAIR  :',Eair_nm1(i_test, j_test), &
1662                    &           ' < ', eair(i_test, j_test), ' < ', Eair_n(i_test, j_test)
1663            ENDIF
1664            WRITE(numout,*) 'TAIR  :',tair_nm1(i_test, j_test), &
1665                 &           ' < ', tair(i_test, j_test), ' < ', tair_n(i_test, j_test)
1666            WRITE(numout,*) 'QAIR  :',qair_nm1(i_test, j_test), &
1667                 &           ' < ', qair(i_test, j_test), ' < ', qair_n(i_test, j_test)
1668            WRITE(numout,*) 'U  :',u_nm1(i_test, j_test), &
1669                 &           ' < ', u(i_test, j_test), ' < ', u_n(i_test, j_test)
1670            WRITE(numout,*) 'V  :',v_nm1(i_test, j_test), &
1671                 &           ' < ', v(i_test, j_test), ' < ', v_n(i_test, j_test)
1672         ENDIF
1673   !---
1674   !--- For precip we suppose that the rain
1675   !--- is the sum over the next 6 hours
1676   !---
1677         IF (itau_split <= nb_spread) THEN
1678            rainf(:,:) = rainf_n(:,:)*(split/REAL(nb_spread))
1679            snowf(:,:) = snowf_n(:,:)*(split/REAL(nb_spread))
1680         ELSE
1681            rainf(:,:) = 0.0
1682            snowf(:,:) = 0.0
1683         ENDIF
1684         IF (check) THEN
1685            WRITE(numout,*) '__ Forcing read at ',itau_split,' :'
1686            WRITE(numout,*) 'Rainf  : ',rainf_nm1(i_test, j_test), &
1687                 &           ' < ', rainf(i_test, j_test), ' < ', rainf_n(i_test, j_test)
1688            WRITE(numout,*) 'Snowf  : ',snowf_nm1(i_test, j_test), &
1689                 &           ' < ', snowf(i_test, j_test), ' < ', snowf_n(i_test, j_test)
1690         ENDIF
1691   !---
1692      ENDIF ! (daily_interpol)
1693   ENDIF
1694!---
1695!--- Here we might put the call to the weather generator ... one day.
1696!--- Pour le moment, le branchement entre interpolation et generateur de temps
1697!--- est fait au-dessus.
1698!---
1699!-   IF ( initialized .AND. weathergen ) THEN
1700!-      ....
1701!-   ENDIF
1702!---
1703!--- At this point the code should be initialized. If not we have a problem !
1704!---
1705   IF ( (itau_read == 0).AND.(itau_split == 0) ) THEN
1706!---
1707      initialized = .TRUE.
1708!---
1709   ELSE
1710      IF ( .NOT. initialized ) THEN
1711         WRITE(numout,*) 'Why is the code forcing_read not initialized ?'
1712         WRITE(numout,*) 'Have you called it with both time-steps set to zero ?'
1713         CALL ipslerr_p(3,'forcing_read_interpol','Pb in initialization','','')
1714      ENDIF
1715   ENDIF
1716
1717END SUBROUTINE forcing_read_interpol
1718
1719
1720!! ==============================================================================================================================\n
1721!! SUBROUTINE   : forcing_just_read
1722!!
1723!>\BRIEF       
1724!!
1725!!\n DESCRIPTION :
1726!!
1727!! RECENT CHANGE(S): None
1728!!
1729!! MAIN OUTPUT VARIABLE(S):
1730!!
1731!! REFERENCE(S) :
1732!!
1733!_ ================================================================================================================================
1734SUBROUTINE forcing_just_read &
1735  & (iim, jjm, zlev, zlev_uv, ttm, itb, ite, &
1736  &  swdown, rainf, snowf, tair, &
1737  &  u, v, qair, pb, lwdown, &
1738  &  SWnet, Eair, petAcoef, peqAcoef, petBcoef, peqBcoef, cdrag, ccanopy, &
1739  &  force_id, wind_N_exists, check)
1740!---------------------------------------------------------------------
1741!- iim        : Size of the grid in x
1742!- jjm        : size of the grid in y
1743!- zlev       : height of the varibales T and Q
1744!- zlev_uv    : height of the varibales U and V
1745!- ttm        : number of time steps in all in the forcing file
1746!- itb, ite   : index of respectively begin and end of read for each variable
1747!- swdown     : Downward solar radiation (W/m^2)
1748!- rainf      : Rainfall (kg/m^2s)
1749!- snowf      : Snowfall (kg/m^2s)
1750!- tair       : 2m air temperature (K)
1751!- u and v    : 2m (in theory !) wind speed (m/s)
1752!- qair       : 2m humidity (kg/kg)
1753!- pb         : Surface pressure (Pa)
1754!- lwdown     : Downward long wave radiation (W/m^2)
1755!-
1756!- From a WATCHOUT file :
1757!- SWnet      : Net surface short-wave flux
1758!- Eair       : Air potential energy
1759!- petAcoef   : Coeficients A from the PBL resolution for T
1760!- peqAcoef   : Coeficients A from the PBL resolution for q
1761!- petBcoef   : Coeficients B from the PBL resolution for T
1762!- peqBcoef   : Coeficients B from the PBL resolution for q
1763!- cdrag      : Surface drag
1764!- ccanopy    : CO2 concentration in the canopy
1765!- force_id   : FLINCOM file id.
1766!-              It is used to close the file at the end of the run.
1767!- wind_N_exists : if Wind_N and Wind_E are in the file (and not just Wind)
1768!- check      : Prompt for reading
1769!---------------------------------------------------------------------
1770   IMPLICIT NONE
1771!-
1772   INTEGER, INTENT(in) :: iim, jjm, ttm
1773   INTEGER, INTENT(in) :: itb, ite
1774   REAL, DIMENSION(iim,jjm), INTENT(out) ::  zlev, zlev_uv, &
1775  &  swdown, rainf, snowf, tair, u, v, qair, pb, lwdown
1776   ! for watchout files
1777   REAL, DIMENSION(iim,jjm), INTENT(out) :: &
1778  &  SWnet, Eair, petAcoef, peqAcoef, petBcoef, peqBcoef, cdrag, ccanopy
1779   INTEGER, INTENT(in) :: force_id
1780!  if Wind_N and Wind_E are in the file (and not just Wind)
1781   LOGICAL, INTENT(in) :: wind_N_exists
1782   LOGICAL :: check
1783   INTEGER :: i, j 
1784   REAL :: rau 
1785
1786!-
1787!---------------------------------------------------------------------
1788   IF ( daily_interpol ) THEN
1789      CALL flinget_buffer (force_id,'Tmin'  , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1790      CALL forcing_zoom(data_full, tair)
1791      CALL flinget_buffer (force_id,'precip' , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1792      CALL forcing_zoom(data_full, rainf)
1793   ELSE
1794      CALL flinget_buffer (force_id,'Tair'  , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1795      CALL forcing_zoom(data_full, tair)
1796      CALL flinget_buffer (force_id,'Snowf' , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1797      CALL forcing_zoom(data_full, snowf)
1798      CALL flinget_buffer (force_id,'Rainf' , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1799      CALL forcing_zoom(data_full, rainf)
1800   ENDIF
1801
1802
1803   CALL flinget_buffer (force_id,'SWdown', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1804   CALL forcing_zoom(data_full, swdown)
1805   CALL flinget_buffer (force_id,'LWdown', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1806   CALL forcing_zoom(data_full, lwdown)
1807
1808   CALL flinget_buffer (force_id,'PSurf' , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1809   CALL forcing_zoom(data_full, pb)
1810   CALL flinget_buffer (force_id,'Qair'  , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1811   CALL forcing_zoom(data_full, qair)
1812!---
1813   IF ( wind_N_exists ) THEN
1814      CALL flinget_buffer (force_id,'Wind_N', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1815      CALL forcing_zoom(data_full, u)
1816      CALL flinget_buffer (force_id,'Wind_E', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1817      CALL forcing_zoom(data_full, v)
1818   ELSE
1819      CALL flinget_buffer (force_id,'Wind',   iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1820      CALL forcing_zoom(data_full, u)
1821      v=0.0
1822   ENDIF
1823
1824!-
1825!- Deal with the height of the atmospheric forcing varibles
1826!-
1827!----
1828   IF ( zheight ) THEN
1829      zlev(:,:) = zlev_fixed 
1830   ELSE IF ( zsigma .OR. zhybrid ) THEN
1831      DO i=1,iim 
1832         DO j=1,jjm 
1833            IF ( tair(i,j) < val_exp ) THEN
1834               rau = pb(i,j)/(cte_molr*tair(i,j)) 
1835                 
1836               zlev(i,j) =  (pb(i,j) - (zhybrid_a + zhybrid_b*pb(i,j)))/(rau * cte_grav) 
1837            ELSE
1838               zlev(i,j) = 0.0 
1839            ENDIF
1840         ENDDO
1841      ENDDO
1842   ELSE IF ( zlevels ) THEN
1843      CALL flinget_buffer (force_id,'Levels', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full) 
1844      CALL forcing_zoom(data_full, zlev) 
1845   ELSE
1846      CALL ipslerr(3, 'forcing_just_read','No case for the vertical levels was specified.', & 
1847           &         'We cannot determine the height for T and Q.','stop readdim2') 
1848   ENDIF
1849   
1850   IF ( zsamelev_uv ) THEN
1851      zlev_uv(:,:) = zlev(:,:) 
1852   ELSE
1853      IF ( zheight ) THEN
1854         zlev_uv(:,:) = zlevuv_fixed 
1855      ELSE IF ( zsigma .OR. zhybrid ) THEN
1856         DO i=1,iim 
1857            DO j=1,jjm 
1858               IF ( tair(i,j) < val_exp ) THEN
1859                  rau = pb(i,j)/(cte_molr*tair(i,j)) 
1860                 
1861                  zlev_uv(i,j) =  (pb(i,j) - (zhybriduv_a + zhybriduv_b*pb(i,j)))/(rau * cte_grav) 
1862               ELSE
1863                  zlev_uv(i,j) = 0.0 
1864               ENDIF
1865            ENDDO
1866         ENDDO
1867      ELSE IF ( zlevels ) THEN
1868         CALL flinget_buffer (force_id,'Levels_uv', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full) 
1869         CALL forcing_zoom(data_full, zlev_uv) 
1870      ELSE
1871         CALL ipslerr(3, 'forcing_just_read','No case for the vertical levels was specified.', & 
1872              &         'We cannot determine the height for U and V.','stop readdim2') 
1873      ENDIF
1874   ENDIF
1875   !----
1876   IF ( is_watchout ) THEN
1877      CALL flinget_buffer (force_id,'levels', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1878      CALL forcing_zoom(data_full, zlev)
1879      !
1880      ! If we are in WATHCOUT it means T,Q are at the same height as U,V
1881      !
1882      zlev_uv(:,:) = zlev(:,:) 
1883      ! Net surface short-wave flux
1884      CALL flinget_buffer (force_id,'SWnet', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1885      CALL forcing_zoom(data_full, SWnet)
1886      ! Air potential energy
1887      CALL flinget_buffer (force_id,'Eair', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1888      CALL forcing_zoom(data_full, Eair)
1889      ! Coeficients A from the PBL resolution for T
1890      CALL flinget_buffer (force_id,'petAcoef', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1891      CALL forcing_zoom(data_full, petAcoef)
1892      ! Coeficients A from the PBL resolution for q
1893      CALL flinget_buffer (force_id,'peqAcoef', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1894      CALL forcing_zoom(data_full, peqAcoef)
1895      ! Coeficients B from the PBL resolution for T
1896      CALL flinget_buffer (force_id,'petBcoef', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1897      CALL forcing_zoom(data_full, petBcoef)
1898      ! Coeficients B from the PBL resolution for q
1899      CALL flinget_buffer (force_id,'peqBcoef', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1900      CALL forcing_zoom(data_full, peqBcoef)
1901      ! Surface drag
1902      CALL flinget_buffer (force_id,'cdrag', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1903      CALL forcing_zoom(data_full, cdrag)
1904      ! CO2 concentration in the canopy
1905      CALL flinget_buffer (force_id,'ccanopy', iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1906      CALL forcing_zoom(data_full, ccanopy)
1907   ENDIF
1908!
1909!----
1910     IF (check) WRITE(numout,*) 'Variables have been extracted between ',itb,' and ',ite,' iterations of the forcing file.'
1911!-------------------------
1912END SUBROUTINE forcing_just_read
1913
1914
1915!! ==============================================================================================================================\n
1916!! SUBROUTINE   : forcing_just_read_tmax
1917!!
1918!>\BRIEF       
1919!!
1920!!\n DESCRIPTION :
1921!!
1922!! RECENT CHANGE(S): None
1923!!
1924!! MAIN OUTPUT VARIABLE(S):
1925!!
1926!! REFERENCE(S) :
1927!!
1928!_ ================================================================================================================================
1929SUBROUTINE forcing_just_read_tmax &
1930  & (iim, jjm, ttm, itb, ite, tmax, force_id )
1931!---------------------------------------------------------------------
1932!- iim        : Size of the grid in x
1933!- jjm        : size of the grid in y
1934!- ttm        : number of time steps in all in the forcing file
1935!- itb, ite   : index of respectively begin and end of read for each variable
1936!- tmax       : 2m air temperature (K)
1937!- force_id   : FLINCOM file id.
1938!-              It is used to close the file at the end of the run.
1939!---------------------------------------------------------------------
1940   IMPLICIT NONE
1941!-
1942   INTEGER, INTENT(in) :: iim, jjm, ttm
1943   INTEGER, INTENT(in) :: itb, ite
1944   REAL, DIMENSION(iim,jjm), INTENT(out) ::  tmax
1945   INTEGER, INTENT(in) :: force_id
1946!-
1947!---------------------------------------------------------------------
1948   CALL flinget_buffer (force_id,'Tmax'  , iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
1949   CALL forcing_zoom(data_full, tmax)
1950!-------------------------
1951END SUBROUTINE forcing_just_read_tmax
1952
1953
1954!! ==============================================================================================================================\n
1955!! SUBROUTINE   : forcing_landind
1956!!
1957!>\BRIEF       
1958!!
1959!!\n DESCRIPTION : This subroutine finds the indices of the land points over which the land
1960!!                 surface scheme is going to run.
1961!!
1962!! RECENT CHANGE(S): None
1963!!
1964!! MAIN OUTPUT VARIABLE(S):
1965!!
1966!! REFERENCE(S) :
1967!!
1968!_ ================================================================================================================================
1969SUBROUTINE forcing_landind(iim, jjm, tair, check, nbindex, kindex, i_test, j_test)
1970!---
1971!---
1972!---
1973  IMPLICIT NONE
1974!-
1975!- ARGUMENTS
1976!-
1977  INTEGER, INTENT(IN) :: iim, jjm
1978  REAL, INTENT(IN)    :: tair(iim,jjm)
1979  INTEGER, INTENT(OUT) :: i_test, j_test, nbindex
1980  INTEGER, INTENT(OUT) :: kindex(iim*jjm)
1981  LOGICAL :: check
1982!-
1983!- LOCAL
1984  INTEGER :: i, j, ig
1985!-
1986!-
1987  ig = 0
1988  i_test = 0
1989  j_test = 0
1990!---
1991  IF (MINVAL(tair(:,:)) < 100.) THEN
1992!----- In this case the 2m temperature is in Celsius
1993     DO j=1,jjm
1994        DO i=1,iim
1995           IF (tair(i,j) < 100.) THEN
1996              ig = ig+1
1997              kindex(ig) = (j-1)*iim+i
1998              !
1999              !  Here we find at random a land-point on which we can do
2000              !  some printouts for test.
2001              !
2002              IF (ig .GT. (iim*jjm)/2 .AND. i_test .LT. 1) THEN
2003                 i_test = i
2004                 j_test = j
2005                 IF (check) THEN
2006                    WRITE(numout,*) 'The test point chosen for output is : ', i_test, j_test
2007                 ENDIF
2008              ENDIF
2009           ENDIF
2010        ENDDO
2011     ENDDO
2012  ELSE 
2013!----- 2m temperature is in Kelvin
2014     DO j=1,jjm
2015        DO i=1,iim
2016           IF (tair(i,j) < 500.) THEN
2017              ig = ig+1
2018              kindex(ig) = (j-1)*iim+i
2019              !
2020              !  Here we find at random a land-point on which we can do
2021              !  some printouts for test.
2022              !
2023              IF (ig .GT. (iim*jjm)/2 .AND. i_test .LT. 1) THEN
2024                 i_test = i
2025                 j_test = j
2026                 IF (check) THEN
2027                    WRITE(numout,*) 'The test point chosen for output is : ', i_test, j_test
2028                 ENDIF
2029              ENDIF
2030           ENDIF
2031        ENDDO
2032     ENDDO
2033  ENDIF
2034
2035  nbindex = ig
2036
2037END SUBROUTINE forcing_landind
2038
2039
2040!! ==============================================================================================================================\n
2041!! SUBROUTINE   : forcing_grid
2042!!
2043!>\BRIEF       
2044!!
2045!!\n DESCRIPTION : This subroutine calculates the longitudes and latitudes of the model grid.
2046!!
2047!! RECENT CHANGE(S): None
2048!!
2049!! MAIN OUTPUT VARIABLE(S):
2050!!
2051!! REFERENCE(S) :
2052!!
2053!_ ================================================================================================================================
2054SUBROUTINE forcing_grid(iim,jjm,llm,lon,lat,init_f)
2055
2056  IMPLICIT NONE
2057
2058  INTEGER, INTENT(in)                   :: iim, jjm, llm
2059  LOGICAL, INTENT(in)                   :: init_f
2060  REAL, DIMENSION(iim,jjm), INTENT(out) :: lon, lat
2061!-
2062  INTEGER :: i,j
2063!-
2064!- Should be unified one day
2065!-
2066  IF ( printlev_loc>=3 ) WRITE(numout,*) 'forcing_grid : options : ', weathergen, interpol
2067!-
2068  IF ( weathergen ) THEN
2069     IF (init_f) THEN
2070        DO i = 1, iim
2071           lon(i,:) = limit_west + merid_res/2. + &
2072                FLOAT(i-1)*(limit_east-limit_west)/FLOAT(iim)
2073        ENDDO
2074        !-
2075        DO j = 1, jjm
2076           lat(:,j) = limit_north - zonal_res/2. - &
2077                FLOAT(j-1)*(limit_north-limit_south)/FLOAT(jjm)
2078        ENDDO
2079     ELSE
2080        IF (is_root_prc) THEN
2081           DO i = 1, iim_g
2082              lon_g(i,:) = limit_west + merid_res/2. + &
2083                   FLOAT(i-1)*(limit_east-limit_west)/FLOAT(iim_g)
2084           ENDDO
2085           !-
2086           DO j = 1, jjm_g
2087              lat_g(:,j) = limit_north - zonal_res/2. - &
2088                   FLOAT(j-1)*(limit_north-limit_south)/FLOAT(jjm_g)
2089           ENDDO
2090        ENDIF
2091        CALL bcast(lon_g)
2092        CALL bcast(lat_g)
2093        lon=lon_g(:,jj_para_begin(mpi_rank):jj_para_end(mpi_rank))
2094        lat=lat_g(:,jj_para_begin(mpi_rank):jj_para_end(mpi_rank))
2095     ENDIF
2096!-
2097  ELSEIF ( interpol ) THEN
2098!-
2099    CALL forcing_zoom(lon_full, lon)
2100    IF ( printlev_loc>=3 ) WRITE(numout,*) 'forcing_grid : out of zoom on lon'
2101    CALL forcing_zoom(lat_full, lat)
2102    IF ( printlev_loc>=3 ) WRITE(numout,*) 'forcing_grid : out of zoom on lat'
2103 
2104 ELSE
2105    CALL ipslerr_p(3,'forcing_grid','Neither interpolation nor weather generator is specified.','','')
2106 ENDIF
2107 
2108 IF ( printlev_loc>=3 ) WRITE(numout,*) 'forcing_grid : ended'
2109 
2110END SUBROUTINE forcing_grid
2111
2112
2113!! ==============================================================================================================================\n
2114!! SUBROUTINE   : forcing_zoom
2115!!
2116!>\BRIEF        This subroutine extracts the region of data over which we wish to run the model.
2117!!
2118!!\n DESCRIPTION : This subroutine extracts the region of data over which we wish to run the model.
2119!!
2120!! RECENT CHANGE(S): None
2121!!
2122!! MAIN OUTPUT VARIABLE(S):
2123!!
2124!! REFERENCE(S) :
2125!!
2126!_ ================================================================================================================================
2127SUBROUTINE forcing_zoom(x_f, x_z)
2128
2129  IMPLICIT NONE
2130!-
2131  REAL, DIMENSION(iim_full, jjm_full), INTENT(IN) :: x_f
2132  REAL, DIMENSION(iim_zoom, jjm_zoom), INTENT(OUT) :: x_z
2133
2134  INTEGER :: i,j
2135!-
2136  DO i=1,iim_zoom
2137     DO j=1,jjm_zoom
2138        x_z(i,j) = x_f(i_index(i),j_index(j))
2139     ENDDO
2140  ENDDO
2141!-
2142END SUBROUTINE forcing_zoom
2143
2144
2145!! ==============================================================================================================================\n
2146!! SUBROUTINE   : forcing_vertical_ioipsl
2147!!
2148!>\BRIEF       
2149!!
2150!!\n DESCRIPTION : This subroutine explores the forcing file to decide what information is available
2151!!                 on the vertical coordinate.
2152!!
2153!! RECENT CHANGE(S): None
2154!!
2155!! MAIN OUTPUT VARIABLE(S):
2156!!
2157!! REFERENCE(S) :
2158!!
2159!_ ================================================================================================================================
2160SUBROUTINE forcing_vertical_ioipsl(force_id)
2161
2162  INTEGER, INTENT(IN) :: force_id
2163
2164  LOGICAL :: var_exists, vara_exists, varb_exists, varuv_exists
2165  LOGICAL :: foundvar
2166  LOGICAL :: levlegacy
2167
2168  !
2169  !- Set all the defaults
2170  !
2171  zfixed=.FALSE.
2172  zsigma=.FALSE.
2173  zhybrid=.FALSE.
2174  zlevels=.FALSE.
2175  zheight=.FALSE.
2176  zsamelev_uv = .TRUE.
2177  levlegacy = .FALSE.
2178  !
2179  foundvar = .FALSE.
2180  !
2181  !- We have a forcing file to explore so let us see if we find any of the conventions
2182  !- which allow us to find the height of T,Q,U and V.
2183  !
2184  IF ( force_id > 0 ) THEN
2185     !
2186     ! Case for sigma levels
2187     !
2188     IF ( .NOT. foundvar ) THEN
2189        CALL flinquery_var(force_id, 'Sigma', var_exists)
2190        CALL flinquery_var(force_id, 'Sigma_uv', varuv_exists)
2191        IF ( var_exists ) THEN
2192           foundvar = .TRUE.
2193           zsigma = .TRUE.
2194           IF ( varuv_exists ) zsamelev_uv = .FALSE.
2195        ENDIF
2196     ENDIF
2197     !
2198     ! Case for Hybrid levels
2199     !
2200     IF ( .NOT. foundvar ) THEN
2201        CALL flinquery_var(force_id, 'HybSigA', vara_exists)
2202        IF ( vara_exists ) THEN
2203           CALL flinquery_var(force_id, 'HybSigB', varb_exists)
2204           IF ( varb_exists ) THEN
2205              var_exists=.TRUE.
2206           ELSE
2207              CALL ipslerr ( 3, 'forcing_vertical_ioipsl','Missing the B coefficient for', &
2208                   &         'Hybrid vertical levels for T and Q','stop readdim2')
2209           ENDIF
2210        ENDIF
2211        CALL flinquery_var(force_id, 'HybSigA_uv', vara_exists)
2212        IF ( vara_exists ) THEN
2213           CALL flinquery_var(force_id, 'HybSigB_uv', varb_exists)
2214           IF ( varb_exists ) THEN
2215              varuv_exists=.TRUE.
2216           ELSE
2217              CALL ipslerr ( 3, 'forcing_vertical_ioipsl','Missing the B coefficient for', &
2218                   &         'Hybrid vertical levels for U and V','stop readdim2')
2219           ENDIF
2220        ENDIF
2221        IF ( var_exists ) THEN
2222           foundvar = .TRUE.
2223           zhybrid = .TRUE.
2224           IF ( varuv_exists ) zsamelev_uv = .FALSE.
2225        ENDIF
2226     ENDIF
2227     !
2228     ! Case for levels (i.e. a 2d time varying field with the height in meters)
2229     !
2230     IF ( .NOT. foundvar ) THEN
2231        CALL flinquery_var(force_id, 'Levels', var_exists)
2232        CALL flinquery_var(force_id, 'Levels_uv', varuv_exists)
2233        IF ( var_exists ) THEN
2234           foundvar = .TRUE.
2235           zlevels = .TRUE.
2236           IF ( varuv_exists ) zsamelev_uv = .FALSE.
2237        ENDIF
2238     ENDIF
2239     !
2240     ! Case where a fixed height is provided in meters
2241     !
2242     IF ( .NOT. foundvar ) THEN
2243        CALL flinquery_var(force_id, 'Height_Lev1', var_exists)
2244        CALL flinquery_var(force_id, 'Height_Levuv', varuv_exists)
2245       IF ( var_exists ) THEN
2246           foundvar = .TRUE.
2247           zheight = .TRUE.
2248           IF ( varuv_exists ) zsamelev_uv = .FALSE.
2249        ENDIF
2250     ENDIF
2251     !
2252     ! Case where a fixed height is provided in meters in the lev variable
2253     !
2254     IF ( .NOT. foundvar ) THEN
2255        CALL flinquery_var(force_id, 'lev', var_exists)
2256        IF ( var_exists ) THEN
2257           foundvar = .TRUE.
2258           zheight = .TRUE.
2259           levlegacy = .TRUE.
2260        ENDIF
2261     ENDIF
2262     !
2263  ENDIF
2264  !
2265  ! We found forcing variables so we need to extract the values if we are dealing with constant values (i.e. all
2266  ! except the case zlevels
2267  !
2268  IF ( foundvar .AND. .NOT. zlevels ) THEN
2269     !
2270     IF ( zheight ) THEN
2271        !
2272        ! Constant height
2273        !
2274        IF ( levlegacy ) THEN
2275           CALL flinget (force_id,'lev',1, 1, 1, 1,  1, 1, zlev_fixed)
2276        ELSE
2277           CALL flinget (force_id,'Height_Lev1',1, 1, 1, 1,  1, 1, zlev_fixed)
2278           IF ( .NOT. zsamelev_uv ) THEN
2279              CALL flinget (force_id,'Height_Levuv',1, 1, 1, 1,  1, 1, zlevuv_fixed)
2280           ENDIF
2281        ENDIF
2282        WRITE(numout,*) "forcing_vertical_ioipsl : case ZLEV : Read from forcing file :", zlev_fixed, zlevuv_fixed
2283        !
2284     ELSE IF ( zsigma .OR. zhybrid ) THEN
2285        !
2286        ! Sigma or hybrid levels
2287        !
2288        IF ( zsigma ) THEN
2289           CALL flinget (force_id,'Sigma',1, 1, 1, 1,  1, 1, zhybrid_b)
2290           zhybrid_a = zero
2291           IF ( .NOT. zsamelev_uv ) THEN
2292              CALL flinget (force_id,'Sigma_uv',1, 1, 1, 1,  1, 1, zhybriduv_b)
2293              zhybriduv_a = zero
2294           ENDIF
2295        ELSE
2296           CALL flinget (force_id,'HybSigB',1, 1, 1, 1,  1, 1, zhybrid_b)
2297           CALL flinget (force_id,'HybSigA',1, 1, 1, 1,  1, 1, zhybrid_a)
2298           IF ( .NOT. zsamelev_uv ) THEN
2299              CALL flinget (force_id,'HybSigB_uv',1, 1, 1, 1,  1, 1, zhybriduv_b)
2300              CALL flinget (force_id,'HybSigA_uv',1, 1, 1, 1,  1, 1, zhybriduv_a)
2301           ENDIF
2302        ENDIF
2303        WRITE(numout,*) "forcing_vertical_ioipsl : case Pressure coordinates : "
2304        WRITE(numout,*) "Read from forcing file :", zhybrid_b, zhybrid_a, zhybriduv_b, zhybriduv_a
2305     ELSE
2306        !
2307        ! Why are we here ???
2308        !
2309        CALL ipslerr ( 3, 'forcing_vertical_ioipsl','What is the option used to describe the height of', &
2310             &         'the atmospheric forcing ?','Please check your forcing file.')
2311     ENDIF
2312  ENDIF
2313  !
2314  !- We have no forcing file to explore or we did not find anything. So revert back to the run.def and
2315  !- read what has been specified by the user.
2316  !
2317  IF ( force_id < 0 .OR. .NOT. foundvar ) THEN
2318     !
2319     !-
2320     !Config Key   = HEIGHT_LEV1
2321     !Config Desc  = Height at which T and Q are given
2322     !Config Def   = 2.0
2323     !Config If    = offline mode
2324     !Config Help  = The atmospheric variables (temperature and specific
2325     !Config         humidity) are measured at a specific level.
2326     !Config         The height of this level is needed to compute
2327     !Config         correctly the turbulent transfer coefficients.
2328     !Config         Look at the description of the forcing
2329     !Config         DATA for the correct value.
2330     !Config Units = [m]
2331     !-
2332     zlev_fixed = 2.0
2333     CALL getin('HEIGHT_LEV1', zlev_fixed)
2334
2335     !-
2336     !Config Key  = HEIGHT_LEVW
2337     !Config Desc = Height at which the wind is given
2338     !Config Def  = 10.0
2339     !Config If   = offline mode
2340     !Config Help = The height at which wind is needed to compute
2341     !Config        correctly the turbulent transfer coefficients.
2342     !Config Units= [m]
2343     !-
2344     zlevuv_fixed = 10.0
2345     CALL getin('HEIGHT_LEVW', zlevuv_fixed)
2346
2347     zheight = .TRUE.
2348
2349     IF ( ABS(zlevuv_fixed-zlev_fixed) > EPSILON(zlev_fixed)) THEN
2350        zsamelev_uv = .FALSE.
2351     ELSE
2352        zsamelev_uv = .TRUE.
2353     ENDIF
2354
2355     CALL ipslerr ( 2, 'forcing_vertical_ioipsl','The height of the atmospheric forcing variables', &
2356          &         'was not found in the netCDF file.','Thus the values in run.def were used ... or their defaults.')
2357  ENDIF
2358
2359END SUBROUTINE forcing_vertical_ioipsl
2360
2361
2362!! ==============================================================================================================================\n
2363!! SUBROUTINE   : domain_size
2364!!
2365!>\BRIEF       
2366!!
2367!!\n DESCRIPTION :
2368!!
2369!! RECENT CHANGE(S): None
2370!!
2371!! MAIN OUTPUT VARIABLE(S):
2372!!
2373!! REFERENCE(S) :
2374!!
2375!_ ================================================================================================================================
2376SUBROUTINE domain_size (limit_west, limit_east, limit_north, limit_south, &
2377     &  iim_f, jjm_f, lon, lat, iim, jjm, iind, jind)
2378 
2379  IMPLICIT NONE
2380  !
2381  ! ARGUMENTS
2382  !
2383  REAL, INTENT(inout)  :: limit_west,limit_east,limit_north,limit_south
2384  INTEGER, INTENT(in)  :: iim_f, jjm_f
2385  REAL, INTENT(in)     :: lon(iim_f, jjm_f), lat(iim_f, jjm_f)
2386  INTEGER, INTENT(out) :: iim,jjm
2387  INTEGER, INTENT(out) :: iind(iim_f), jind(jjm_f)
2388  !
2389  ! LOCAL
2390  !
2391  INTEGER :: i, j
2392  REAL :: lolo
2393  LOGICAL :: over_dateline = .FALSE.
2394  !
2395  !
2396  IF ( ( ABS(limit_east) .GT. 180. ) .OR. &
2397       ( ABS(limit_west) .GT. 180. ) ) THEN
2398     WRITE(numout,*) 'Limites Ouest, Est: ',limit_west,limit_east
2399     CALL ipslerr_p (3,'domain_size', &
2400 &        'Longitudes problem.','In run.def file :', &
2401 &        'limit_east > 180. or limit_west > 180.')
2402  ENDIF
2403  !
2404  IF ( limit_west .GT. limit_east ) over_dateline = .TRUE.
2405  !
2406  IF ( ( limit_south .LT. -90. ) .OR. &
2407       ( limit_north .GT. 90. ) .OR. &
2408       ( limit_south .GE. limit_north ) ) THEN
2409     WRITE(numout,*) 'Limites Nord, Sud: ',limit_north,limit_south
2410     CALL ipslerr_p (3,'domain_size', &
2411 &        'Latitudes problem.','In run.def file :', &
2412 &        'limit_south < -90. or limit_north > 90. or limit_south >= limit_north')
2413  ENDIF
2414  !
2415  ! Here we assume that the grid of the forcing data is regular. Else we would have
2416  ! to do more work to find the index table.
2417  !
2418  iim = 0
2419  DO i=1,iim_f
2420     !
2421     lolo =  lon(i,1)
2422     IF ( lon(i,1) .GT. 180. ) lolo =  lon(i,1) - 360.
2423     IF ( lon(i,1) .LT. -180. ) lolo =  lon(i,1) + 360.
2424     !
2425     IF (lon(i,1) < limit_west) iim_g_begin = i+1
2426     IF (lon(i,1) < limit_east) iim_g_end = i
2427     !
2428     IF ( over_dateline ) THEN
2429        IF ( lolo .LE. limit_west .OR. lolo .GE. limit_east ) THEN
2430           iim = iim + 1
2431           iind(iim) = i
2432        ENDIF
2433     ELSE
2434        IF ( lolo .GE. limit_west .AND. lolo .LE. limit_east ) THEN
2435           iim = iim + 1
2436           iind(iim) = i
2437        ENDIF
2438     ENDIF
2439     !
2440  ENDDO
2441  !
2442  jjm = 0
2443  DO j=1,jjm_f
2444     IF (lat(1,j) > limit_north) jjm_g_begin = j+1
2445     IF (lat(1,j) > limit_south) jjm_g_end = j
2446     !
2447     IF ( lat(1,j) .GE. limit_south .AND. lat(1,j) .LE. limit_north) THEN
2448        jjm = jjm + 1
2449        jind(jjm) = j
2450     ENDIF
2451  ENDDO
2452
2453  WRITE(numout,*) 'Domain zoom size: iim, jjm = ', iim, jjm
2454
2455  END SUBROUTINE domain_size
2456
2457
2458!! ==============================================================================================================================\n
2459!! SUBROUTINE   : flinget_buffer
2460!!
2461!>\BRIEF       
2462!!
2463!!\n DESCRIPTION :  This subroutine is a wrap of flinget/IOIPSL. The arguments are the same.
2464!!                  flinget_buffer will call flinget and buffer the forcing data localy in this subroutine.
2465!!                  According to the variable NBUFF set in run.def, several time steps can be read at the same time
2466!!                  from the forcing file. If NBUFF=0, the full forcing file is read.
2467!!                  The output, data_full, from this subroutine is always only one time step of corresponding to itb.
2468!!                  itb must be equal to ite.
2469!!
2470!! RECENT CHANGE(S): None
2471!!
2472!! MAIN OUTPUT VARIABLE(S):
2473!!
2474!! REFERENCE(S) :
2475!!
2476!_ ================================================================================================================================
2477  SUBROUTINE flinget_buffer(force_id, varname, iim_full, jjm_full, llm_full, ttm, itb, ite, data_full)
2478   
2479    !! Input arguments
2480    INTEGER, INTENT(in)          :: force_id                      !! Id for forcing file
2481    CHARACTER(len=*), INTENT(in) :: varname                       !! Name of current variable to be read
2482    INTEGER, INTENT(in)          :: iim_full, jjm_full, llm_full  !! Horizontal and vertical domaine
2483    INTEGER, INTENT(in)          :: ttm                           !! Full lenght of forcing file
2484    INTEGER, INTENT(in)          :: itb, ite                      !! Time step to be read from forcing file. itb must be equal to ite
2485
2486    !! Output argument
2487    REAL, DIMENSION(iim_full, jjm_full), INTENT(out) :: data_full !! Data for time step itb.
2488
2489    !! Define specific type to buffer data together with name and index
2490    TYPE buffer_type
2491       CHARACTER(len=20) :: name                   !! Name of variable in forcing file
2492       INTEGER :: istart                           !! Start index of current buffered data
2493       INTEGER :: iend                             !! End index of current buffered data
2494       REAL, ALLOCATABLE, DIMENSION(:,:,:) :: data !! Data read from forcing file for intervall [istart,iend]
2495    END TYPE buffer_type
2496   
2497    !! Local variables
2498    INTEGER, PARAMETER :: maxvar=20                          !! Max number of variables to be buffered
2499    TYPE(buffer_type), DIMENSION(maxvar),SAVE :: data_buffer !! Containing all variables and the current buffered data
2500    INTEGER, SAVE   :: nbuff                                 !! Number of time steps to be buffered
2501    INTEGER, SAVE   :: lastindex=0                           !! Current number of variables stored in data_buffer
2502    INTEGER, SAVE   :: ttm0                                  !! Time lenght of forcing file, stored for test purpose
2503    LOGICAL, SAVE   :: first=.TRUE.                          !! First call to this subroutine
2504    INTEGER         :: index                                 !! Index in data_buffer for current variable
2505    INTEGER         :: i, ierr                               !! Loop and error variables
2506   
2507   
2508    !! 1. Initialization
2509    IF (first) THEN
2510       data_buffer(:)%name='undef'
2511       ! Read NBUFF from run.def
2512       ! Note that getin_p is not used because this subroutine might be called only by master process
2513
2514       !Config Key  = NBUFF
2515       !Config Desc = Number of time steps of data to buffer between each reading of the forcing file
2516       !Config If   = OFF_LINE
2517       !Config Help = The full simulation time length will be read if NBUFF equal 0.
2518       !Config        NBUFF > 1 can be used for smaller regions or site simulations only.
2519       !Config Def  = 1
2520       !Config Units= -
2521
2522       nbuff=1
2523       CALL getin('NBUFF', nbuff)
2524
2525       IF (nbuff == 0 .OR. nbuff >ttm) THEN
2526          ! Set nbuff as the full forcing file lenght
2527          nbuff=ttm
2528       ELSE IF (nbuff < 0) THEN
2529          ! Negativ nbuff not possible
2530          CALL ipslerr_p(3,'flinget_buffer','NBUFF must be a positiv number','Set NBUFF=0 for full simulation lenght','')
2531       END IF
2532       WRITE(numout,*)'flinget_buffer: NBUFF=',nbuff,' number of time step will be buffered'
2533       WRITE(numout,*)'flinget_buffer: Choose a lower value for NBUFF if problem with memory'
2534       
2535       ! Save dimensions to check following timesteps
2536       ! ttm is the full lenght of forcing file
2537       ttm0=ttm
2538       
2539       first=.FALSE.
2540    END IF
2541
2542    !! 2. Coeherence tests on input arguments
2543    IF (ttm /= ttm0) THEN
2544       WRITE(numout,*)'Problem with ttm=',ttm,' ttm0=',ttm0
2545       CALL ipslerr_p(3,'flinget_buffer','Error with ttm and ttm0','','')
2546    END IF
2547    IF (itb /= ite) THEN
2548       WRITE(numout,*) 'There is a problem. Why is itb not equal ite ?'
2549       WRITE(numout,*) 'itb=',itb,' ite=',ite,' varname=',varname
2550       CALL ipslerr_p(3,'flinget_buffer','ite not equal itb','','')
2551    END IF
2552
2553
2554    !! 3. Find index for current variable
2555    index=0
2556    DO i=1, maxvar
2557       IF ( trim(varname) == data_buffer(i)%name ) THEN
2558          index=i
2559          CYCLE
2560       END IF
2561    END DO
2562   
2563    !! 4. Initialize and allocate if necesary the current variable
2564    IF ( index == 0 ) THEN
2565       ! The variable was not found
2566       ! This must be the first time for current variable
2567       index=lastindex+1
2568       lastindex=index
2569       IF (index > maxvar) CALL ipslerr_p(3,'flinget_buffer','to many variables','maxvar is too small','')
2570       
2571       ! Initialize the data_buffer for this index
2572       data_buffer(index)%name=trim(varname)
2573       ALLOCATE(data_buffer(index)%data(iim_full,jjm_full,nbuff),stat=ierr)
2574       IF (ierr /= 0) CALL ipslerr_p(3,'flinget_buffer','pb alloc data_buffer%data','for variable=',varname)
2575       data_buffer(index)%istart=0
2576       data_buffer(index)%iend=0
2577    END IF
2578   
2579   
2580    !! 5. Call flinget if current time step (itb) is outside the buffered intervall
2581    IF (( itb > data_buffer(index)%iend ) .OR. ( itb < data_buffer(index)%istart )) THEN
2582       ! itb is not in the time slice previously read or it is the first time to read
2583       ! Reading of forcing file will now be done
2584       ! First recalculate index to be read
2585       data_buffer(index)%istart = itb
2586       data_buffer(index)%iend   = itb + nbuff - 1
2587       
2588!       WRITE(numout,*) 'Now do flinget for ',varname,', itb=',itb,', istart=',&
2589!            data_buffer(index)%istart,', iend=',data_buffer(index)%iend
2590       CALL flinget (force_id,varname, iim_full, jjm_full, llm_full, ttm, data_buffer(index)%istart, &
2591            data_buffer(index)%iend, data_buffer(index)%data(:,:,:))
2592    END IF
2593   
2594    !! 6. Initialize the output variable with data from buffered variable
2595    ! Find index for the time step corrsponding to itb in the time slice previously read from forcing file
2596    i=itb-data_buffer(index)%istart+1
2597    ! Initialize output variable
2598    data_full(:,:) = data_buffer(index)%data(:,:,i)
2599   
2600   
2601  END SUBROUTINE flinget_buffer
2602
2603END MODULE readdim2
Note: See TracBrowser for help on using the repository browser.