source: LMDZ6/trunk/libf/phylmd/dyn1d/mod_1D_cases_read2.f90 @ 5327

Last change on this file since 5327 was 5302, checked in by abarral, 2 weeks ago

Turn compar1d.h date_cas.h into module
Move fcg_racmo.h to obsolete

File size: 60.4 KB
Line 
1!
2! $Id: mod_1D_cases_read.F90 2373 2015-10-13 17:28:01Z jyg $
3!
4MODULE mod_1D_cases_read2
5  USE netcdf, ONLY: nf90_get_var, nf90_noerr, nf90_inq_varid, nf90_inquire_dimension, nf90_strerror, nf90_open, &
6          nf90_nowrite, nf90_inq_dimid
7!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
8  !Declarations specifiques au cas standard
9  character*80 :: fich_cas
10  ! Discr?tisation
11  integer nlev_cas, nt_cas
12
13
14  !profils environnementaux
15  real, allocatable::  plev_cas(:,:),plevh_cas(:)
16  real, allocatable::  ap_cas(:),bp_cas(:)
17
18  real, allocatable::  z_cas(:,:),zh_cas(:)
19  real, allocatable::  t_cas(:,:),q_cas(:,:),qv_cas(:,:),ql_cas(:,:),qi_cas(:,:),rh_cas(:,:)
20  real, allocatable::  th_cas(:,:),thv_cas(:,:),thl_cas(:,:),rv_cas(:,:)
21  real, allocatable::  u_cas(:,:),v_cas(:,:),vitw_cas(:,:),omega_cas(:,:)
22
23  !forcing
24  real, allocatable::  ht_cas(:,:),vt_cas(:,:),dt_cas(:,:),dtrad_cas(:,:)
25  real, allocatable::  hth_cas(:,:),vth_cas(:,:),dth_cas(:,:)
26  real, allocatable::  hq_cas(:,:),vq_cas(:,:),dq_cas(:,:)
27  real, allocatable::  hr_cas(:,:),vr_cas(:,:),dr_cas(:,:)
28  real, allocatable::  hu_cas(:,:),vu_cas(:,:),du_cas(:,:)
29  real, allocatable::  hv_cas(:,:),vv_cas(:,:),dv_cas(:,:)
30  real, allocatable::  ug_cas(:,:),vg_cas(:,:)
31  real, allocatable::  lat_cas(:),sens_cas(:),ts_cas(:),ps_cas(:),ustar_cas(:)
32  real, allocatable::  uw_cas(:,:),vw_cas(:,:),q1_cas(:,:),q2_cas(:,:),tke_cas(:)
33
34  !champs interpoles
35  real, allocatable::  plev_prof_cas(:)
36  real, allocatable::  t_prof_cas(:)
37  real, allocatable::  theta_prof_cas(:)
38  real, allocatable::  thl_prof_cas(:)
39  real, allocatable::  thv_prof_cas(:)
40  real, allocatable::  q_prof_cas(:)
41  real, allocatable::  qv_prof_cas(:)
42  real, allocatable::  ql_prof_cas(:)
43  real, allocatable::  qi_prof_cas(:)
44  real, allocatable::  rh_prof_cas(:)
45  real, allocatable::  rv_prof_cas(:)
46  real, allocatable::  u_prof_cas(:)
47  real, allocatable::  v_prof_cas(:)       
48  real, allocatable::  vitw_prof_cas(:)
49  real, allocatable::  omega_prof_cas(:)
50  real, allocatable::  ug_prof_cas(:)
51  real, allocatable::  vg_prof_cas(:)
52  real, allocatable::  ht_prof_cas(:)
53  real, allocatable::  hth_prof_cas(:)
54  real, allocatable::  hq_prof_cas(:)
55  real, allocatable::  vt_prof_cas(:)
56  real, allocatable::  vth_prof_cas(:)
57  real, allocatable::  vq_prof_cas(:)
58  real, allocatable::  dt_prof_cas(:)
59  real, allocatable::  dth_prof_cas(:)
60  real, allocatable::  dtrad_prof_cas(:)
61  real, allocatable::  dq_prof_cas(:)
62  real, allocatable::  hu_prof_cas(:)
63  real, allocatable::  hv_prof_cas(:)
64  real, allocatable::  vu_prof_cas(:)
65  real, allocatable::  vv_prof_cas(:)
66  real, allocatable::  du_prof_cas(:)
67  real, allocatable::  dv_prof_cas(:)
68  real, allocatable::  uw_prof_cas(:)
69  real, allocatable::  vw_prof_cas(:)
70  real, allocatable::  q1_prof_cas(:)
71  real, allocatable::  q2_prof_cas(:)
72
73
74  real lat_prof_cas,sens_prof_cas,ts_prof_cas,ps_prof_cas,ustar_prof_cas,tke_prof_cas
75  real o3_cas,orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,heat_rough,rugos_cas,sand_cas,clay_cas
76
77
78
79CONTAINS
80
81  SUBROUTINE read_1D_cas
82    implicit none
83
84    INTEGER nid,rid,ierr
85    INTEGER ii,jj
86
87    fich_cas='setup/cas.nc'
88    print*,'fich_cas ',fich_cas
89    ierr = nf90_open(fich_cas,nf90_nowrite,nid)
90    print*,'fich_cas,nf90_nowrite,nid ',fich_cas,nf90_nowrite,nid
91    if (ierr.NE.nf90_noerr) then
92       write(*,*) 'ERROR: GROS Pb opening forcings nc file '
93       write(*,*) nf90_strerror(ierr)
94       stop ""
95    endif
96    !.......................................................................
97    ierr=nf90_inq_dimid(nid,'lat',rid)
98    IF (ierr.NE.nf90_noerr) THEN
99       print*, 'Oh probleme lecture dimension lat'
100    ENDIF
101    ierr=nf90_inquire_dimension(nid,rid,len=ii)
102    print*,'OK1 nid,rid,lat',nid,rid,ii
103    !.......................................................................
104    ierr=nf90_inq_dimid(nid,'lon',rid)
105    IF (ierr.NE.nf90_noerr) THEN
106       print*, 'Oh probleme lecture dimension lon'
107    ENDIF
108    ierr=nf90_inquire_dimension(nid,rid,len=jj)
109    print*,'OK2 nid,rid,lat',nid,rid,jj
110    !.......................................................................
111    ierr=nf90_inq_dimid(nid,'lev',rid)
112    IF (ierr.NE.nf90_noerr) THEN
113       print*, 'Oh probleme lecture dimension zz'
114    ENDIF
115    ierr=nf90_inquire_dimension(nid,rid,len=nlev_cas)
116    print*,'OK3 nid,rid,nlev_cas',nid,rid,nlev_cas
117    !.......................................................................
118    ierr=nf90_inq_dimid(nid,'time',rid)
119    print*,'nid,rid',nid,rid
120    nt_cas=0
121    IF (ierr.NE.nf90_noerr) THEN
122       stop 'probleme lecture dimension sens'
123    ENDIF
124    ierr=nf90_inquire_dimension(nid,rid,len=nt_cas)
125    print*,'OK4 nid,rid,nt_cas',nid,rid,nt_cas
126
127!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
128    !profils moyens:
129    allocate(plev_cas(nlev_cas,nt_cas))       
130    allocate(z_cas(nlev_cas,nt_cas))
131    allocate(t_cas(nlev_cas,nt_cas),q_cas(nlev_cas,nt_cas),rh_cas(nlev_cas,nt_cas))
132    allocate(th_cas(nlev_cas,nt_cas),rv_cas(nlev_cas,nt_cas))
133    allocate(u_cas(nlev_cas,nt_cas))
134    allocate(v_cas(nlev_cas,nt_cas))
135
136    !forcing
137    allocate(ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas),dt_cas(nlev_cas,nt_cas),dtrad_cas(nlev_cas,nt_cas))
138    allocate(hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas),dq_cas(nlev_cas,nt_cas))
139    allocate(hth_cas(nlev_cas,nt_cas),vth_cas(nlev_cas,nt_cas),dth_cas(nlev_cas,nt_cas))
140    allocate(hr_cas(nlev_cas,nt_cas),vr_cas(nlev_cas,nt_cas),dr_cas(nlev_cas,nt_cas))
141    allocate(hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas),du_cas(nlev_cas,nt_cas))
142    allocate(hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas),dv_cas(nlev_cas,nt_cas))
143    allocate(vitw_cas(nlev_cas,nt_cas))
144    allocate(ug_cas(nlev_cas,nt_cas))
145    allocate(vg_cas(nlev_cas,nt_cas))
146    allocate(lat_cas(nt_cas),sens_cas(nt_cas),ts_cas(nt_cas),ps_cas(nt_cas),ustar_cas(nt_cas))
147    allocate(uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas),q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas))
148
149
150    !champs interpoles
151    allocate(plev_prof_cas(nlev_cas))
152    allocate(t_prof_cas(nlev_cas))
153    allocate(q_prof_cas(nlev_cas))
154    allocate(u_prof_cas(nlev_cas))
155    allocate(v_prof_cas(nlev_cas))
156
157    allocate(vitw_prof_cas(nlev_cas))
158    allocate(ug_prof_cas(nlev_cas))
159    allocate(vg_prof_cas(nlev_cas))
160    allocate(ht_prof_cas(nlev_cas))
161    allocate(hq_prof_cas(nlev_cas))
162    allocate(hu_prof_cas(nlev_cas))
163    allocate(hv_prof_cas(nlev_cas))
164    allocate(vt_prof_cas(nlev_cas))
165    allocate(vq_prof_cas(nlev_cas))
166    allocate(vu_prof_cas(nlev_cas))
167    allocate(vv_prof_cas(nlev_cas))
168    allocate(dt_prof_cas(nlev_cas))
169    allocate(dtrad_prof_cas(nlev_cas))
170    allocate(dq_prof_cas(nlev_cas))
171    allocate(du_prof_cas(nlev_cas))
172    allocate(dv_prof_cas(nlev_cas))
173    allocate(uw_prof_cas(nlev_cas))
174    allocate(vw_prof_cas(nlev_cas))
175    allocate(q1_prof_cas(nlev_cas))
176    allocate(q2_prof_cas(nlev_cas))
177
178    print*,'Allocations OK'
179    call read_cas2(nid,nlev_cas,nt_cas                                       &
180         ,z_cas,plev_cas,t_cas,q_cas,rh_cas,th_cas,rv_cas,u_cas,v_cas         &
181         ,ug_cas,vg_cas,vitw_cas,du_cas,hu_cas,vu_cas,dv_cas,hv_cas,vv_cas    &
182         ,dt_cas,dtrad_cas,ht_cas,vt_cas,dq_cas,hq_cas,vq_cas                 &
183         ,dth_cas,hth_cas,vth_cas,dr_cas,hr_cas,vr_cas,sens_cas,lat_cas,ts_cas&
184         ,ustar_cas,uw_cas,vw_cas,q1_cas,q2_cas)
185    print*,'Read cas OK'
186
187
188  END SUBROUTINE read_1D_cas
189  !**********************************************************************************************
190  SUBROUTINE read2_1D_cas
191    implicit none
192
193    INTEGER nid,rid,ierr
194    INTEGER ii,jj
195
196    fich_cas='setup/cas.nc'
197    print*,'fich_cas ',fich_cas
198    ierr = nf90_open(fich_cas,nf90_nowrite,nid)
199    print*,'fich_cas,nf90_nowrite,nid ',fich_cas,nf90_nowrite,nid
200    if (ierr.NE.nf90_noerr) then
201       write(*,*) 'ERROR: GROS Pb opening forcings nc file '
202       write(*,*) nf90_strerror(ierr)
203       stop ""
204    endif
205    !.......................................................................
206    ierr=nf90_inq_dimid(nid,'lat',rid)
207    IF (ierr.NE.nf90_noerr) THEN
208       print*, 'Oh probleme lecture dimension lat'
209    ENDIF
210    ierr=nf90_inquire_dimension(nid,rid,len=ii)
211    print*,'OK1 read2: nid,rid,lat',nid,rid,ii
212    !.......................................................................
213    ierr=nf90_inq_dimid(nid,'lon',rid)
214    IF (ierr.NE.nf90_noerr) THEN
215       print*, 'Oh probleme lecture dimension lon'
216    ENDIF
217    ierr=nf90_inquire_dimension(nid,rid,len=jj)
218    print*,'OK2 read2: nid,rid,lat',nid,rid,jj
219    !.......................................................................
220    ierr=nf90_inq_dimid(nid,'nlev',rid)
221    IF (ierr.NE.nf90_noerr) THEN
222       print*, 'Oh probleme lecture dimension nlev'
223    ENDIF
224    ierr=nf90_inquire_dimension(nid,rid,len=nlev_cas)
225    print*,'OK3 read2: nid,rid,nlev_cas',nid,rid,nlev_cas
226    !.......................................................................
227    ierr=nf90_inq_dimid(nid,'time',rid)
228    nt_cas=0
229    IF (ierr.NE.nf90_noerr) THEN
230       stop 'Oh probleme lecture dimension time'
231    ENDIF
232    ierr=nf90_inquire_dimension(nid,rid,len=nt_cas)
233    print*,'OK4 read2: nid,rid,nt_cas',nid,rid,nt_cas
234
235!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
236    !profils moyens:
237    allocate(plev_cas(nlev_cas,nt_cas),plevh_cas(nlev_cas+1))       
238    allocate(z_cas(nlev_cas,nt_cas),zh_cas(nlev_cas+1))
239    allocate(ap_cas(nlev_cas+1),bp_cas(nlev_cas+1))
240    allocate(t_cas(nlev_cas,nt_cas),q_cas(nlev_cas,nt_cas),qv_cas(nlev_cas,nt_cas),ql_cas(nlev_cas,nt_cas), &
241         qi_cas(nlev_cas,nt_cas),rh_cas(nlev_cas,nt_cas))
242    allocate(th_cas(nlev_cas,nt_cas),thl_cas(nlev_cas,nt_cas),thv_cas(nlev_cas,nt_cas),rv_cas(nlev_cas,nt_cas))
243    allocate(u_cas(nlev_cas,nt_cas),v_cas(nlev_cas,nt_cas),vitw_cas(nlev_cas,nt_cas),omega_cas(nlev_cas,nt_cas))
244
245    !forcing
246    allocate(ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas),dt_cas(nlev_cas,nt_cas),dtrad_cas(nlev_cas,nt_cas))
247    allocate(hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas),dq_cas(nlev_cas,nt_cas))
248    allocate(hth_cas(nlev_cas,nt_cas),vth_cas(nlev_cas,nt_cas),dth_cas(nlev_cas,nt_cas))
249    allocate(hr_cas(nlev_cas,nt_cas),vr_cas(nlev_cas,nt_cas),dr_cas(nlev_cas,nt_cas))
250    allocate(hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas),du_cas(nlev_cas,nt_cas))
251    allocate(hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas),dv_cas(nlev_cas,nt_cas))
252    allocate(ug_cas(nlev_cas,nt_cas))
253    allocate(vg_cas(nlev_cas,nt_cas))
254    allocate(lat_cas(nt_cas),sens_cas(nt_cas),ts_cas(nt_cas),ps_cas(nt_cas),ustar_cas(nt_cas),tke_cas(nt_cas))
255    allocate(uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas),q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas))
256
257
258
259    !champs interpoles
260    allocate(plev_prof_cas(nlev_cas))
261    allocate(t_prof_cas(nlev_cas))
262    allocate(theta_prof_cas(nlev_cas))
263    allocate(thl_prof_cas(nlev_cas))
264    allocate(thv_prof_cas(nlev_cas))
265    allocate(q_prof_cas(nlev_cas))
266    allocate(qv_prof_cas(nlev_cas))
267    allocate(ql_prof_cas(nlev_cas))
268    allocate(qi_prof_cas(nlev_cas))
269    allocate(rh_prof_cas(nlev_cas))
270    allocate(rv_prof_cas(nlev_cas))
271    allocate(u_prof_cas(nlev_cas))
272    allocate(v_prof_cas(nlev_cas))
273    allocate(vitw_prof_cas(nlev_cas))
274    allocate(omega_prof_cas(nlev_cas))
275    allocate(ug_prof_cas(nlev_cas))
276    allocate(vg_prof_cas(nlev_cas))
277    allocate(ht_prof_cas(nlev_cas))
278    allocate(hth_prof_cas(nlev_cas))
279    allocate(hq_prof_cas(nlev_cas))
280    allocate(hu_prof_cas(nlev_cas))
281    allocate(hv_prof_cas(nlev_cas))
282    allocate(vt_prof_cas(nlev_cas))
283    allocate(vth_prof_cas(nlev_cas))
284    allocate(vq_prof_cas(nlev_cas))
285    allocate(vu_prof_cas(nlev_cas))
286    allocate(vv_prof_cas(nlev_cas))
287    allocate(dt_prof_cas(nlev_cas))
288    allocate(dth_prof_cas(nlev_cas))
289    allocate(dtrad_prof_cas(nlev_cas))
290    allocate(dq_prof_cas(nlev_cas))
291    allocate(du_prof_cas(nlev_cas))
292    allocate(dv_prof_cas(nlev_cas))
293    allocate(uw_prof_cas(nlev_cas))
294    allocate(vw_prof_cas(nlev_cas))
295    allocate(q1_prof_cas(nlev_cas))
296    allocate(q2_prof_cas(nlev_cas))
297
298    print*,'Allocations OK'
299    call read2_cas (nid,nlev_cas,nt_cas,                                                                     &
300         ap_cas,bp_cas,z_cas,plev_cas,zh_cas,plevh_cas,t_cas,th_cas,thv_cas,thl_cas,qv_cas,                    &
301         ql_cas,qi_cas,rh_cas,rv_cas,u_cas,v_cas,vitw_cas,omega_cas,ug_cas,vg_cas,du_cas,hu_cas,vu_cas,        &
302         dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dq_cas,hq_cas,vq_cas,dth_cas,hth_cas,vth_cas,               &
303         dr_cas,hr_cas,vr_cas,dtrad_cas,sens_cas,lat_cas,ts_cas,ps_cas,ustar_cas,tke_cas,                      &
304         uw_cas,vw_cas,q1_cas,q2_cas,orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,heat_rough, &
305         o3_cas,rugos_cas,clay_cas,sand_cas)
306    print*,'Read2 cas OK'
307    do ii=1,nlev_cas
308       print*,'apres read2_cas, plev_cas=',ii,plev_cas(ii,1)
309    enddo
310
311
312  END SUBROUTINE read2_1D_cas
313
314  !**********************************************************************************************
315  SUBROUTINE old_read_SCM_cas
316    use netcdf, only: nf90_get_var
317    USE date_cas_mod_h
318    implicit none
319
320    INTEGER nid,rid,ierr
321    INTEGER ii,jj,timeid
322    REAL, ALLOCATABLE :: time_val(:)
323
324    fich_cas='cas.nc'
325    print*,'fich_cas ',fich_cas
326    ierr = nf90_open(fich_cas,nf90_nowrite,nid)
327    print*,'fich_cas,nf90_nowrite,nid ',fich_cas,nf90_nowrite,nid
328    if (ierr.NE.nf90_noerr) then
329       write(*,*) 'ERROR: GROS Pb opening forcings nc file '
330       write(*,*) nf90_strerror(ierr)
331       stop ""
332    endif
333    !.......................................................................
334    ierr=nf90_inq_dimid(nid,'lat',rid)
335    IF (ierr.NE.nf90_noerr) THEN
336       print*, 'Oh probleme lecture dimension lat'
337    ENDIF
338    ierr=nf90_inquire_dimension(nid,rid,len=ii)
339    print*,'OK1 read2: nid,rid,lat',nid,rid,ii
340    !.......................................................................
341    ierr=nf90_inq_dimid(nid,'lon',rid)
342    IF (ierr.NE.nf90_noerr) THEN
343       print*, 'Oh probleme lecture dimension lon'
344    ENDIF
345    ierr=nf90_inquire_dimension(nid,rid,len=jj)
346    print*,'OK2 read2: nid,rid,lat',nid,rid,jj
347    !.......................................................................
348    ierr=nf90_inq_dimid(nid,'lev',rid)
349    IF (ierr.NE.nf90_noerr) THEN
350       print*, 'Oh probleme lecture dimension nlev'
351    ENDIF
352    ierr=nf90_inquire_dimension(nid,rid,len=nlev_cas)
353    print*,'OK3 read2: nid,rid,nlev_cas',nid,rid,nlev_cas
354    IF ( .NOT. ( nlev_cas > 10 .AND. nlev_cas < 1000 )) THEN
355       print*,'Valeur de nlev_cas peu probable'
356       STOP
357    ENDIF
358    !.......................................................................
359    ierr=nf90_inq_dimid(nid,'time',rid)
360    nt_cas=0
361    IF (ierr.NE.nf90_noerr) THEN
362       stop 'Oh probleme lecture dimension time'
363    ENDIF
364    ierr=nf90_inquire_dimension(nid,rid,len=nt_cas)
365    print*,'OK4 read2: nid,rid,nt_cas',nid,rid,nt_cas
366    ! Lecture de l'axe des temps
367    print*,'LECTURE DU TEMPS'
368    ierr=nf90_inq_varid(nid,'time',timeid)
369    if(ierr/=nf90_noerr) then
370       print *,'Variable time manquante dans cas.nc:'
371       ierr=nf90_noerr
372    else
373       allocate(time_val(nt_cas))
374       ierr = NF90_GET_VAR(nid,timeid,time_val)
375       if(ierr/=nf90_noerr) then
376          print *,'Pb a la lecture de time cas.nc: '
377       endif
378    endif
379    IF (nt_cas>1) THEN
380       pdt_cas=time_val(2)-time_val(1)
381    ELSE
382       pdt_cas=0.
383    ENDIF
384
385
386!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
387    !profils moyens:
388    allocate(plev_cas(nlev_cas,nt_cas),plevh_cas(nlev_cas+1))       
389    allocate(z_cas(nlev_cas,nt_cas),zh_cas(nlev_cas+1))
390    allocate(ap_cas(nlev_cas+1),bp_cas(nlev_cas+1))
391    allocate(t_cas(nlev_cas,nt_cas),q_cas(nlev_cas,nt_cas),qv_cas(nlev_cas,nt_cas),ql_cas(nlev_cas,nt_cas), &
392         qi_cas(nlev_cas,nt_cas),rh_cas(nlev_cas,nt_cas))
393    allocate(th_cas(nlev_cas,nt_cas),thl_cas(nlev_cas,nt_cas),thv_cas(nlev_cas,nt_cas),rv_cas(nlev_cas,nt_cas))
394    allocate(u_cas(nlev_cas,nt_cas),v_cas(nlev_cas,nt_cas),vitw_cas(nlev_cas,nt_cas),omega_cas(nlev_cas,nt_cas))
395
396    !forcing
397    allocate(ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas),dt_cas(nlev_cas,nt_cas),dtrad_cas(nlev_cas,nt_cas))
398    allocate(hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas),dq_cas(nlev_cas,nt_cas))
399    allocate(hth_cas(nlev_cas,nt_cas),vth_cas(nlev_cas,nt_cas),dth_cas(nlev_cas,nt_cas))
400    allocate(hr_cas(nlev_cas,nt_cas),vr_cas(nlev_cas,nt_cas),dr_cas(nlev_cas,nt_cas))
401    allocate(hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas),du_cas(nlev_cas,nt_cas))
402    allocate(hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas),dv_cas(nlev_cas,nt_cas))
403    allocate(ug_cas(nlev_cas,nt_cas))
404    allocate(vg_cas(nlev_cas,nt_cas))
405    allocate(lat_cas(nt_cas),sens_cas(nt_cas),ts_cas(nt_cas),ps_cas(nt_cas),ustar_cas(nt_cas),tke_cas(nt_cas))
406    allocate(uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas),q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas))
407
408
409
410    !champs interpoles
411    allocate(plev_prof_cas(nlev_cas))
412    allocate(t_prof_cas(nlev_cas))
413    allocate(theta_prof_cas(nlev_cas))
414    allocate(thl_prof_cas(nlev_cas))
415    allocate(thv_prof_cas(nlev_cas))
416    allocate(q_prof_cas(nlev_cas))
417    allocate(qv_prof_cas(nlev_cas))
418    allocate(ql_prof_cas(nlev_cas))
419    allocate(qi_prof_cas(nlev_cas))
420    allocate(rh_prof_cas(nlev_cas))
421    allocate(rv_prof_cas(nlev_cas))
422    allocate(u_prof_cas(nlev_cas))
423    allocate(v_prof_cas(nlev_cas))
424    allocate(vitw_prof_cas(nlev_cas))
425    allocate(omega_prof_cas(nlev_cas))
426    allocate(ug_prof_cas(nlev_cas))
427    allocate(vg_prof_cas(nlev_cas))
428    allocate(ht_prof_cas(nlev_cas))
429    allocate(hth_prof_cas(nlev_cas))
430    allocate(hq_prof_cas(nlev_cas))
431    allocate(hu_prof_cas(nlev_cas))
432    allocate(hv_prof_cas(nlev_cas))
433    allocate(vt_prof_cas(nlev_cas))
434    allocate(vth_prof_cas(nlev_cas))
435    allocate(vq_prof_cas(nlev_cas))
436    allocate(vu_prof_cas(nlev_cas))
437    allocate(vv_prof_cas(nlev_cas))
438    allocate(dt_prof_cas(nlev_cas))
439    allocate(dth_prof_cas(nlev_cas))
440    allocate(dtrad_prof_cas(nlev_cas))
441    allocate(dq_prof_cas(nlev_cas))
442    allocate(du_prof_cas(nlev_cas))
443    allocate(dv_prof_cas(nlev_cas))
444    allocate(uw_prof_cas(nlev_cas))
445    allocate(vw_prof_cas(nlev_cas))
446    allocate(q1_prof_cas(nlev_cas))
447    allocate(q2_prof_cas(nlev_cas))
448
449    print*,'Allocations OK'
450    call old_read_SCM (nid,nlev_cas,nt_cas,                                                                     &
451         ap_cas,bp_cas,z_cas,plev_cas,zh_cas,plevh_cas,t_cas,th_cas,thv_cas,thl_cas,qv_cas,                    &
452         ql_cas,qi_cas,rh_cas,rv_cas,u_cas,v_cas,vitw_cas,omega_cas,ug_cas,vg_cas,du_cas,hu_cas,vu_cas,        &
453         dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dq_cas,hq_cas,vq_cas,dth_cas,hth_cas,vth_cas,               &
454         dr_cas,hr_cas,vr_cas,dtrad_cas,sens_cas,lat_cas,ts_cas,ps_cas,ustar_cas,tke_cas,                      &
455         uw_cas,vw_cas,q1_cas,q2_cas,orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,heat_rough, &
456         o3_cas,rugos_cas,clay_cas,sand_cas)
457    print*,'Read2 cas OK'
458    do ii=1,nlev_cas
459       print*,'apres read2_cas, plev_cas=',ii,plev_cas(ii,1)
460    enddo
461
462
463  END SUBROUTINE old_read_SCM_cas
464
465
466!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
467  SUBROUTINE deallocate2_1D_cases
468    !profils environnementaux:
469    deallocate(plev_cas,plevh_cas)
470
471    deallocate(z_cas,zh_cas)
472    deallocate(ap_cas,bp_cas)
473    deallocate(t_cas,q_cas,qv_cas,ql_cas,qi_cas,rh_cas)
474    deallocate(th_cas,thl_cas,thv_cas,rv_cas)
475    deallocate(u_cas,v_cas,vitw_cas,omega_cas)
476
477    !forcing
478    deallocate(ht_cas,vt_cas,dt_cas,dtrad_cas)
479    deallocate(hq_cas,vq_cas,dq_cas)
480    deallocate(hth_cas,vth_cas,dth_cas)
481    deallocate(hr_cas,vr_cas,dr_cas)
482    deallocate(hu_cas,vu_cas,du_cas)
483    deallocate(hv_cas,vv_cas,dv_cas)
484    deallocate(ug_cas)
485    deallocate(vg_cas)
486    deallocate(lat_cas,sens_cas,ts_cas,ps_cas,ustar_cas,tke_cas,uw_cas,vw_cas,q1_cas,q2_cas)
487
488    !champs interpoles
489    deallocate(plev_prof_cas)
490    deallocate(t_prof_cas)
491    deallocate(theta_prof_cas)
492    deallocate(thl_prof_cas)
493    deallocate(thv_prof_cas)
494    deallocate(q_prof_cas)
495    deallocate(qv_prof_cas)
496    deallocate(ql_prof_cas)
497    deallocate(qi_prof_cas)
498    deallocate(rh_prof_cas)
499    deallocate(rv_prof_cas)
500    deallocate(u_prof_cas)
501    deallocate(v_prof_cas)
502    deallocate(vitw_prof_cas)
503    deallocate(omega_prof_cas)
504    deallocate(ug_prof_cas)
505    deallocate(vg_prof_cas)
506    deallocate(ht_prof_cas)
507    deallocate(hq_prof_cas)
508    deallocate(hu_prof_cas)
509    deallocate(hv_prof_cas)
510    deallocate(vt_prof_cas)
511    deallocate(vq_prof_cas)
512    deallocate(vu_prof_cas)
513    deallocate(vv_prof_cas)
514    deallocate(dt_prof_cas)
515    deallocate(dtrad_prof_cas)
516    deallocate(dq_prof_cas)
517    deallocate(du_prof_cas)
518    deallocate(dv_prof_cas)
519    deallocate(t_prof_cas)
520    deallocate(u_prof_cas)
521    deallocate(v_prof_cas)
522    deallocate(uw_prof_cas)
523    deallocate(vw_prof_cas)
524    deallocate(q1_prof_cas)
525    deallocate(q2_prof_cas)
526
527  END SUBROUTINE deallocate2_1D_cases
528
529
530END MODULE mod_1D_cases_read2
531!=====================================================================
532subroutine read_cas2(nid,nlevel,ntime                          &
533     ,zz,pp,temp,qv,rh,theta,rv,u,v,ug,vg,w,                   &
534     du,hu,vu,dv,hv,vv,dt,dtrad,ht,vt,dq,hq,vq,                     &
535     dth,hth,vth,dr,hr,vr,sens,flat,ts,ustar,uw,vw,q1,q2)
536
537  !program reading forcing of the case study
538  USE netcdf, ONLY: nf90_get_var, nf90_noerr, nf90_inq_varid, nf90_inquire_dimension, nf90_strerror, nf90_open, &
539          nf90_nowrite, nf90_inq_dimid
540  implicit none
541
542  integer ntime,nlevel
543
544  real zz(nlevel,ntime)
545  real pp(nlevel,ntime)
546  real temp(nlevel,ntime),qv(nlevel,ntime),rh(nlevel,ntime)
547  real theta(nlevel,ntime),rv(nlevel,ntime)
548  real u(nlevel,ntime)
549  real v(nlevel,ntime)
550  real ug(nlevel,ntime)
551  real vg(nlevel,ntime)
552  real w(nlevel,ntime)
553  real du(nlevel,ntime),hu(nlevel,ntime),vu(nlevel,ntime)
554  real dv(nlevel,ntime),hv(nlevel,ntime),vv(nlevel,ntime)
555  real dt(nlevel,ntime),ht(nlevel,ntime),vt(nlevel,ntime)
556  real dtrad(nlevel,ntime)
557  real dq(nlevel,ntime),hq(nlevel,ntime),vq(nlevel,ntime)
558  real dth(nlevel,ntime),hth(nlevel,ntime),vth(nlevel,ntime)
559  real dr(nlevel,ntime),hr(nlevel,ntime),vr(nlevel,ntime)
560  real flat(ntime),sens(ntime),ts(ntime),ustar(ntime)
561  real uw(nlevel,ntime),vw(nlevel,ntime),q1(nlevel,ntime),q2(nlevel,ntime),resul(nlevel,ntime),resul1(ntime)
562
563
564  integer nid, ierr, ierr1,ierr2,rid,i
565  integer nbvar3d
566  parameter(nbvar3d=39)
567  integer var3didin(nbvar3d)
568  character*5 name_var(1:nbvar3d)
569  data name_var/'zz','pp','temp','qv','rh','theta','rv','u','v','ug','vg','w','advu','hu','vu',&
570       'advv','hv','vv','advT','hT','vT','advq','hq','vq','advth','hth','vth','advr','hr','vr',&
571       'radT','uw','vw','q1','q2','sens','flat','ts','ustar'/
572
573
574  do i=1,nbvar3d
575     print *,'Dans read_cas2, on va lire ',nid,i,name_var(i)
576  enddo
577  do i=1,nbvar3d
578     ierr=nf90_inq_varid(nid,name_var(i),var3didin(i))
579     print *,'ierr=',i,ierr,name_var(i),var3didin(i)
580     if(ierr/=nf90_noerr) then
581        print *,'Variable manquante dans cas.nc:',name_var(i)
582     endif
583  enddo
584  do i=1,nbvar3d
585     print *,'Dans read_cas2, on va lire ',var3didin(i),name_var(i)
586     if(i.LE.35) then
587        ierr = NF90_GET_VAR(nid,var3didin(i),resul, count = [1, 1, nlevel, ntime])
588        print *,'Dans read_cas2, on a lu ',ierr,var3didin(i),name_var(i)
589        if(ierr/=nf90_noerr) then
590           print *,'Pb a la lecture de cas.nc: ',name_var(i)
591           stop "getvarup"
592        endif
593     else
594        print *,'Dans read_cas2, on a lu ',ierr,var3didin(i),name_var(i)
595        ierr = NF90_GET_VAR(nid,var3didin(i),resul1, count = [1, 1, ntime])
596        if(ierr/=nf90_noerr) then
597           print *,'Pb a la lecture de cas.nc: ',name_var(i)
598           stop "getvarup"
599        endif
600     endif
601     select case(i)
602     case(1) ; zz=resul
603     case(2) ; pp=resul
604     case(3) ; temp=resul
605     case(4) ; qv=resul
606     case(5) ; rh=resul
607     case(6) ; theta=resul
608     case(7) ; rv=resul
609     case(8) ; u=resul
610     case(9) ; v=resul
611     case(10) ; ug=resul
612     case(11) ; vg=resul
613     case(12) ; w=resul
614     case(13) ; du=resul
615     case(14) ; hu=resul
616     case(15) ; vu=resul
617     case(16) ; dv=resul
618     case(17) ; hv=resul
619     case(18) ; vv=resul
620     case(19) ; dt=resul
621     case(20) ; ht=resul
622     case(21) ; vt=resul
623     case(22) ; dq=resul
624     case(23) ; hq=resul
625     case(24) ; vq=resul
626     case(25) ; dth=resul
627     case(26) ; hth=resul
628     case(27) ; vth=resul
629     case(28) ; dr=resul
630     case(29) ; hr=resul
631     case(30) ; vr=resul
632     case(31) ; dtrad=resul
633     case(32) ; uw=resul
634     case(33) ; vw=resul
635     case(34) ; q1=resul
636     case(35) ; q2=resul
637     case(36) ; sens=resul1
638     case(37) ; flat=resul1
639     case(38) ; ts=resul1
640     case(39) ; ustar=resul1
641     end select
642  enddo
643
644  return
645end subroutine read_cas2
646!======================================================================
647subroutine read2_cas(nid,nlevel,ntime,                                       &
648     ap,bp,zz,pp,zzh,pph,temp,theta,thv,thl,qv,ql,qi,rh,rv,u,v,vitw,omega,ug,vg,&
649     du,hu,vu,dv,hv,vv,dt,ht,vt,dq,hq,vq,                                    &
650     dth,hth,vth,dr,hr,vr,dtrad,sens,flat,ts,ps,ustar,tke,uw,vw,q1,q2,       &
651     orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,          &
652     heat_rough,o3_cas,rugos_cas,clay_cas,sand_cas)
653
654  !program reading forcing of the case study
655  USE netcdf, ONLY: nf90_get_var, nf90_noerr, nf90_inq_varid, nf90_inquire_dimension, nf90_strerror, nf90_open, &
656          nf90_nowrite, nf90_inq_dimid
657  implicit none
658
659  integer ntime,nlevel
660
661  real ap(nlevel+1),bp(nlevel+1)
662  real zz(nlevel,ntime),zzh(nlevel+1)
663  real pp(nlevel,ntime),pph(nlevel+1)
664  real temp(nlevel,ntime),qv(nlevel,ntime),ql(nlevel,ntime),qi(nlevel,ntime),rh(nlevel,ntime)
665  real theta(nlevel,ntime),thv(nlevel,ntime),thl(nlevel,ntime),rv(nlevel,ntime)
666  real u(nlevel,ntime),v(nlevel,ntime)
667  real ug(nlevel,ntime),vg(nlevel,ntime)
668  real vitw(nlevel,ntime),omega(nlevel,ntime)
669  real du(nlevel,ntime),hu(nlevel,ntime),vu(nlevel,ntime)
670  real dv(nlevel,ntime),hv(nlevel,ntime),vv(nlevel,ntime)
671  real dt(nlevel,ntime),ht(nlevel,ntime),vt(nlevel,ntime)
672  real dtrad(nlevel,ntime)
673  real dq(nlevel,ntime),hq(nlevel,ntime),vq(nlevel,ntime)
674  real dth(nlevel,ntime),hth(nlevel,ntime),vth(nlevel,ntime),hthl(nlevel,ntime)
675  real dr(nlevel,ntime),hr(nlevel,ntime),vr(nlevel,ntime)
676  real flat(ntime),sens(ntime),ustar(ntime)
677  real uw(nlevel,ntime),vw(nlevel,ntime),q1(nlevel,ntime),q2(nlevel,ntime)
678  real ts(ntime),ps(ntime),tke(ntime)
679  real orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,heat_rough,o3_cas,rugos_cas,clay_cas,sand_cas
680  real apbp(nlevel+1),resul(nlevel,ntime),resul1(nlevel),resul2(ntime),resul3
681
682
683  integer nid, ierr,ierr1,ierr2,rid,i
684  integer nbvar3d
685  parameter(nbvar3d=62)
686  integer var3didin(nbvar3d),missing_var(nbvar3d)
687  character*12 name_var(1:nbvar3d)
688  data name_var/'coor_par_a','coor_par_b','height_h','pressure_h',&
689       'w','omega','ug','vg','uadv','uadvh','uadvv','vadv','vadvh','vadvv','tadv','tadvh','tadvv',&
690       'qadv','qadvh','qadvv','thadv','thadvh','thadvv','thladvh','radv','radvh','radvv','radcool','q1','q2','ustress','vstress', &
691       'rh',&
692       'height_f','pressure_f','temp','theta','thv','thl','qv','ql','qi','rv','u','v',&
693       'sfc_sens_flx','sfc_lat_flx','ts','ps','ustar','tke',&
694       'orog','albedo','emiss','t_skin','q_skin','mom_rough','heat_rough','o3','rugos','clay','sand'/
695  do i=1,nbvar3d
696     missing_var(i)=0.
697  enddo
698
699  !-----------------------------------------------------------------------
700  do i=1,nbvar3d
701     ierr=nf90_inq_varid(nid,name_var(i),var3didin(i))
702     if(ierr/=nf90_noerr) then
703        print *,'Variable manquante dans cas.nc:',i,name_var(i)
704        ierr=nf90_noerr
705        missing_var(i)=1
706     else
707        !-----------------------------------------------------------------------
708        if(i.LE.4) then     ! Lecture des coord pression en (nlevelp1,lat,lon)
709           ierr = NF90_GET_VAR(nid,var3didin(i),apbp, count = [1, 1, nlevel + 1])
710           print *,'read2_cas(apbp), on a lu ',i,name_var(i)
711           if(ierr/=nf90_noerr) then
712              print *,'Pb a la lecture de cas.nc: ',name_var(i)
713              stop "getvarup"
714           endif
715           !-----------------------------------------------------------------------
716        else if(i.gt.4.and.i.LE.45) then   ! Lecture des variables en (time,nlevel,lat,lon)
717           ierr = NF90_GET_VAR(nid,var3didin(i),resul, count = [1, 1, nlevel, ntime])
718           print *,'read2_cas(resul), on a lu ',i,name_var(i)
719           if(ierr/=nf90_noerr) then
720              print *,'Pb a la lecture de cas.nc: ',name_var(i)
721              stop "getvarup"
722           endif
723           !-----------------------------------------------------------------------
724        else if (i.gt.45.and.i.LE.51) then   ! Lecture des variables en (time,lat,lon)
725           ierr = NF90_GET_VAR(nid,var3didin(i),resul2, count = [1, 1, ntime])
726           print *,'read2_cas(resul2), on a lu ',i,name_var(i)
727           if(ierr/=nf90_noerr) then
728              print *,'Pb a la lecture de cas.nc: ',name_var(i)
729              stop "getvarup"
730           endif
731           !-----------------------------------------------------------------------
732        else     ! Lecture des constantes (lat,lon)
733           ierr = NF90_GET_VAR(nid,var3didin(i),resul3)
734           print *,'read2_cas(resul3), on a lu ',i,name_var(i)
735           if(ierr/=nf90_noerr) then
736              print *,'Pb a la lecture de cas.nc: ',name_var(i)
737              stop "getvarup"
738           endif
739        endif
740     endif
741     !-----------------------------------------------------------------------
742     select case(i)
743     case(1) ; ap=apbp       ! donnees indexees en nlevel+1
744     case(2) ; bp=apbp
745     case(3) ; zzh=apbp
746     case(4) ; pph=apbp
747     case(5) ; vitw=resul    ! donnees indexees en nlevel,time
748     case(6) ; omega=resul
749     case(7) ; ug=resul
750     case(8) ; vg=resul
751     case(9) ; du=resul
752     case(10) ; hu=resul
753     case(11) ; vu=resul
754     case(12) ; dv=resul
755     case(13) ; hv=resul
756     case(14) ; vv=resul
757     case(15) ; dt=resul
758     case(16) ; ht=resul
759     case(17) ; vt=resul
760     case(18) ; dq=resul
761     case(19) ; hq=resul
762     case(20) ; vq=resul
763     case(21) ; dth=resul
764     case(22) ; hth=resul
765     case(23) ; vth=resul
766     case(24) ; hthl=resul
767     case(25) ; dr=resul
768     case(26) ; hr=resul
769     case(27) ; vr=resul
770     case(28) ; dtrad=resul
771     case(29) ; q1=resul
772     case(30) ; q2=resul
773     case(31) ; uw=resul
774     case(32) ; vw=resul
775     case(33) ; rh=resul
776     case(34) ; zz=resul      ! donnees en time,nlevel pour profil initial
777     case(35) ; pp=resul
778     case(36) ; temp=resul
779     case(37) ; theta=resul
780     case(38) ; thv=resul
781     case(39) ; thl=resul
782     case(40) ; qv=resul
783     case(41) ; ql=resul
784     case(42) ; qi=resul
785     case(43) ; rv=resul
786     case(44) ; u=resul
787     case(45) ; v=resul
788     case(46) ; sens=resul2   ! donnees indexees en time
789     case(47) ; flat=resul2
790     case(48) ; ts=resul2
791     case(49) ; ps=resul2
792     case(50) ; ustar=resul2
793     case(51) ; tke=resul2
794     case(52) ; orog_cas=resul3      ! constantes
795     case(53) ; albedo_cas=resul3
796     case(54) ; emiss_cas=resul3
797     case(55) ; t_skin_cas=resul3
798     case(56) ; q_skin_cas=resul3
799     case(57) ; mom_rough=resul3
800     case(58) ; heat_rough=resul3
801     case(59) ; o3_cas=resul3       
802     case(60) ; rugos_cas=resul3
803     case(61) ; clay_cas=resul3
804     case(62) ; sand_cas=resul3
805     end select
806     resul=0.
807     resul1=0.
808     resul2=0.
809     resul3=0.
810  enddo
811  !-----------------------------------------------------------------------
812
813
814  return
815end subroutine read2_cas
816
817!======================================================================
818subroutine old_read_SCM(nid,nlevel,ntime,                                       &
819     ap,bp,zz,pp,zzh,pph,temp,theta,thv,thl,qv,ql,qi,rh,rv,u,v,vitw,omega,ug,vg,&
820     du,hu,vu,dv,hv,vv,dt,ht,vt,dq,hq,vq,                                    &
821     dth,hth,vth,dr,hr,vr,dtrad,sens,flat,ts,ps,ustar,tke,uw,vw,q1,q2,       &
822     orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,          &
823     heat_rough,o3_cas,rugos_cas,clay_cas,sand_cas)
824
825  !program reading forcing of the case study
826  USE netcdf, ONLY: nf90_get_var, nf90_noerr, nf90_inq_varid, nf90_inquire_dimension, nf90_strerror, nf90_open, &
827          nf90_nowrite, nf90_inq_dimid
828  implicit none
829
830  integer ntime,nlevel,k,t
831
832  real ap(nlevel+1),bp(nlevel+1)
833  real zz(nlevel,ntime),zzh(nlevel+1)
834  real pp(nlevel,ntime),pph(nlevel+1)
835  !profils initiaux
836  real temp0(nlevel),qv0(nlevel),ql0(nlevel),qi0(nlevel),u0(nlevel),v0(nlevel),tke0(nlevel)
837  real pp0(nlevel)   
838  real temp(nlevel,ntime),qv(nlevel,ntime),ql(nlevel,ntime),qi(nlevel,ntime),rh(nlevel,ntime)
839  real theta(nlevel,ntime),thv(nlevel,ntime),thl(nlevel,ntime),rv(nlevel,ntime)
840  real u(nlevel,ntime),v(nlevel,ntime),tke(nlevel,ntime)
841  real ug(nlevel,ntime),vg(nlevel,ntime)
842  real vitw(nlevel,ntime),omega(nlevel,ntime)
843  real du(nlevel,ntime),hu(nlevel,ntime),vu(nlevel,ntime)
844  real dv(nlevel,ntime),hv(nlevel,ntime),vv(nlevel,ntime)
845  real dt(nlevel,ntime),ht(nlevel,ntime),vt(nlevel,ntime)
846  real dtrad(nlevel,ntime)
847  real dq(nlevel,ntime),hq(nlevel,ntime),vq(nlevel,ntime)
848  real dth(nlevel,ntime),hth(nlevel,ntime),vth(nlevel,ntime),hthl(nlevel,ntime)
849  real dr(nlevel,ntime),hr(nlevel,ntime),vr(nlevel,ntime)
850  real flat(ntime),sens(ntime),ustar(ntime)
851  real uw(nlevel,ntime),vw(nlevel,ntime),q1(nlevel,ntime),q2(nlevel,ntime)
852  real ts(ntime),ps(ntime)
853  real orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,heat_rough,o3_cas,rugos_cas,clay_cas,sand_cas
854  real apbp(nlevel+1),resul(nlevel,ntime),resul1(nlevel),resul2(ntime),resul3
855
856
857  integer nid, ierr,ierr1,ierr2,rid,i
858  integer nbvar3d
859  parameter(nbvar3d=70)
860  integer var3didin(nbvar3d),missing_var(nbvar3d)
861  character*13 name_var(1:nbvar3d)
862  data name_var/'coor_par_a','coor_par_b','height_h','pressure_h',&
863       'temp','qv','ql','qi','u','v','tke','pressure',&
864       'w','omega','ug','vg','uadv','uadvh','uadvv','vadv','vadvh','vadvv','tadv','tadvh','tadvv',&
865       'qvadv','qvadvh','qvadvv','thadv','thadvh','thadvv','thladvh','radv','radvh','radvv','radcool','q1','q2','ustress', &
866       'vstress','rh',&
867       'height_f','pressure_forc','tempt','theta','thv','thl','qvt','qlt','qit','rv','ut','vt','tket',&
868       'sfc_sens_flx','sfc_lat_flx','ts','ps','ustar',&
869       'orog','albedo','emiss','t_skin','q_skin','mom_rough','heat_rough','o3','rugos','clay','sand'/
870  do i=1,nbvar3d
871     missing_var(i)=0.
872  enddo
873
874  !-----------------------------------------------------------------------
875
876  print*,'ON EST LA'
877  do i=1,nbvar3d
878     ierr=nf90_inq_varid(nid,name_var(i),var3didin(i))
879     if(ierr/=nf90_noerr) then
880        print *,'Variable manquante dans cas.nc:',i,name_var(i)
881        ierr=nf90_noerr
882        missing_var(i)=1
883     else
884        !-----------------------------------------------------------------------
885        if(i.LE.4) then     ! Lecture des coord pression en (nlevelp1,lat,lon)
886           ierr = NF90_GET_VAR(nid,var3didin(i),apbp)
887           print *,'read2_cas(apbp), on a lu ',i,name_var(i)
888           if(ierr/=nf90_noerr) then
889              print *,'Pb a la lecture de cas.nc: ',name_var(i)
890              stop "getvarup"
891           endif
892           !-----------------------------------------------------------------------
893        else if(i.gt.4.and.i.LE.12) then   ! Lecture des variables en (time,nlevel,lat,lon)
894           ierr = NF90_GET_VAR(nid,var3didin(i),resul1)
895           print *,'read2_cas(resul1), on a lu ',i,name_var(i)
896           if(ierr/=nf90_noerr) then
897              print *,'Pb a la lecture de cas.nc: ',name_var(i)
898              stop "getvarup"
899           endif
900           print*,'Lecture de la variable #i ',i,name_var(i),minval(resul1),maxval(resul1)
901           !-----------------------------------------------------------------------
902        else if(i.gt.12.and.i.LE.54) then   ! Lecture des variables en (time,nlevel,lat,lon)
903           ierr = NF90_GET_VAR(nid,var3didin(i),resul)
904           print *,'read2_cas(resul), on a lu ',i,name_var(i)
905           if(ierr/=nf90_noerr) then
906              print *,'Pb a la lecture de cas.nc: ',name_var(i)
907              stop "getvarup"
908           endif
909           print*,'Lecture de la variable #i ',i,name_var(i),minval(resul),maxval(resul)
910           !-----------------------------------------------------------------------
911        else if (i.gt.54.and.i.LE.65) then   ! Lecture des variables en (time,lat,lon)
912           ierr = NF90_GET_VAR(nid,var3didin(i),resul2)
913           print *,'read2_cas(resul2), on a lu ',i,name_var(i)
914           if(ierr/=nf90_noerr) then
915              print *,'Pb a la lecture de cas.nc: ',name_var(i)
916              stop "getvarup"
917           endif
918           print*,'Lecture de la variable #i  ',i,name_var(i),minval(resul2),maxval(resul2)
919           !-----------------------------------------------------------------------
920        else     ! Lecture des constantes (lat,lon)
921           ierr = NF90_GET_VAR(nid,var3didin(i),resul3)
922           print *,'read2_cas(resul3), on a lu ',i,name_var(i)
923           if(ierr/=nf90_noerr) then
924              print *,'Pb a la lecture de cas.nc: ',name_var(i)
925              stop "getvarup"
926           endif
927           print*,'Lecture de la variable #i ',i,name_var(i),resul3
928        endif
929     endif
930     !-----------------------------------------------------------------------
931     select case(i)
932        !case(1) ; ap=apbp       ! donnees indexees en nlevel+1
933        ! case(2) ; bp=apbp
934     case(3) ; zzh=apbp
935     case(4) ; pph=apbp
936     case(5) ; temp0=resul1    ! donnees initiales
937     case(6) ; qv0=resul1
938     case(7) ; ql0=resul1
939     case(8) ; qi0=resul1
940     case(9) ; u0=resul1
941     case(10) ; v0=resul1
942     case(11) ; tke0=resul1
943     case(12) ; pp0=resul1
944     case(13) ; vitw=resul    ! donnees indexees en nlevel,time
945     case(14) ; omega=resul
946     case(15) ; ug=resul
947     case(16) ; vg=resul
948     case(17) ; du=resul
949     case(18) ; hu=resul
950     case(19) ; vu=resul
951     case(20) ; dv=resul
952     case(21) ; hv=resul
953     case(22) ; vv=resul
954     case(23) ; dt=resul
955     case(24) ; ht=resul
956     case(25) ; vt=resul
957     case(26) ; dq=resul
958     case(27) ; hq=resul
959     case(28) ; vq=resul
960     case(29) ; dth=resul
961     case(30) ; hth=resul
962     case(31) ; vth=resul
963     case(32) ; hthl=resul
964     case(33) ; dr=resul
965     case(34) ; hr=resul
966     case(35) ; vr=resul
967     case(36) ; dtrad=resul
968     case(37) ; q1=resul
969     case(38) ; q2=resul
970     case(39) ; uw=resul
971     case(40) ; vw=resul
972     case(41) ; rh=resul
973     case(42) ; zz=resul      ! donnees en time,nlevel pour profil initial
974     case(43) ; pp=resul
975     case(44) ; temp=resul
976     case(45) ; theta=resul
977     case(46) ; thv=resul
978     case(47) ; thl=resul
979     case(48) ; qv=resul
980     case(49) ; ql=resul
981     case(50) ; qi=resul
982     case(51) ; rv=resul
983     case(52) ; u=resul
984     case(53) ; v=resul
985     case(54) ; tke=resul
986     case(55) ; sens=resul2   ! donnees indexees en time
987     case(56) ; flat=resul2
988     case(57) ; ts=resul2
989     case(58) ; ps=resul2
990     case(59) ; ustar=resul2
991     case(60) ; orog_cas=resul3      ! constantes
992     case(61) ; albedo_cas=resul3
993     case(62) ; emiss_cas=resul3
994     case(63) ; t_skin_cas=resul3
995     case(64) ; q_skin_cas=resul3
996     case(65) ; mom_rough=resul3
997     case(66) ; heat_rough=resul3
998     case(67) ; o3_cas=resul3       
999     case(68) ; rugos_cas=resul3
1000     case(69) ; clay_cas=resul3
1001     case(70) ; sand_cas=resul3
1002     end select
1003     resul=0.
1004     resul1=0.
1005     resul2=0.
1006     resul3=0.
1007  enddo
1008  print*,'Lecture de la variable APRES ,sens ',minval(sens),maxval(sens)
1009  print*,'Lecture de la variable APRES ,flat ',minval(flat),maxval(flat)
1010
1011  !CR:ATTENTION EN ATTENTE DE REGLER LA QUESTION DU PAS DE TEMPS INITIAL
1012  do t=1,ntime
1013     do k=1,nlevel
1014        temp(k,t)=temp0(k)
1015        qv(k,t)=qv0(k)
1016        ql(k,t)=ql0(k)
1017        qi(k,t)=qi0(k)
1018        u(k,t)=u0(k)
1019        v(k,t)=v0(k)
1020        tke(k,t)=tke0(k)
1021     enddo
1022  enddo
1023  !-----------------------------------------------------------------------
1024
1025  return
1026end subroutine old_read_SCM
1027!======================================================================
1028
1029!======================================================================
1030SUBROUTINE interp_case_time2(day,day1,annee_ref                &
1031     !    &         ,year_cas,day_cas,nt_cas,pdt_forc,nlev_cas      &
1032     ,nt_cas,nlev_cas                                       &
1033     ,ts_cas,ps_cas,plev_cas,t_cas,q_cas,u_cas,v_cas               &
1034     ,ug_cas,vg_cas,vitw_cas,du_cas,hu_cas,vu_cas           &
1035     ,dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dtrad_cas   &
1036     ,dq_cas,hq_cas,vq_cas,lat_cas,sens_cas,ustar_cas       &
1037     ,uw_cas,vw_cas,q1_cas,q2_cas                           &
1038     ,ts_prof_cas,plev_prof_cas,t_prof_cas,q_prof_cas       &
1039     ,u_prof_cas,v_prof_cas,ug_prof_cas,vg_prof_cas         &
1040     ,vitw_prof_cas,du_prof_cas,hu_prof_cas,vu_prof_cas     &
1041     ,dv_prof_cas,hv_prof_cas,vv_prof_cas,dt_prof_cas       &
1042     ,ht_prof_cas,vt_prof_cas,dtrad_prof_cas,dq_prof_cas    &
1043     ,hq_prof_cas,vq_prof_cas,lat_prof_cas,sens_prof_cas    &
1044     ,ustar_prof_cas,uw_prof_cas,vw_prof_cas,q1_prof_cas,q2_prof_cas)
1045  USE compar1d_mod_h
1046  USE date_cas_mod_h
1047  implicit none
1048
1049  !---------------------------------------------------------------------------------------
1050  ! Time interpolation of a 2D field to the timestep corresponding to day
1051  !
1052  ! day: current julian day (e.g. 717538.2)
1053  ! day1: first day of the simulation
1054  ! nt_cas: total nb of data in the forcing
1055  ! pdt_cas: total time interval (in sec) between 2 forcing data
1056  !---------------------------------------------------------------------------------------
1057
1058  ! inputs:
1059  integer annee_ref
1060  integer nt_cas,nlev_cas
1061  real day, day1,day_cas
1062  real ts_cas(nt_cas),ps_cas(nt_cas)
1063  real plev_cas(nlev_cas,nt_cas)
1064  real t_cas(nlev_cas,nt_cas),q_cas(nlev_cas,nt_cas)
1065  real u_cas(nlev_cas,nt_cas),v_cas(nlev_cas,nt_cas)
1066  real ug_cas(nlev_cas,nt_cas),vg_cas(nlev_cas,nt_cas)
1067  real vitw_cas(nlev_cas,nt_cas)
1068  real du_cas(nlev_cas,nt_cas),hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas)
1069  real dv_cas(nlev_cas,nt_cas),hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas)
1070  real dt_cas(nlev_cas,nt_cas),ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas)
1071  real dtrad_cas(nlev_cas,nt_cas)
1072  real dq_cas(nlev_cas,nt_cas),hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas)
1073  real lat_cas(nt_cas)
1074  real sens_cas(nt_cas)
1075  real ustar_cas(nt_cas),uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas)
1076  real q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas)
1077
1078  ! outputs:
1079  real plev_prof_cas(nlev_cas)
1080  real t_prof_cas(nlev_cas),q_prof_cas(nlev_cas)
1081  real u_prof_cas(nlev_cas),v_prof_cas(nlev_cas)
1082  real ug_prof_cas(nlev_cas),vg_prof_cas(nlev_cas)
1083  real vitw_prof_cas(nlev_cas)
1084  real du_prof_cas(nlev_cas),hu_prof_cas(nlev_cas),vu_prof_cas(nlev_cas)
1085  real dv_prof_cas(nlev_cas),hv_prof_cas(nlev_cas),vv_prof_cas(nlev_cas)
1086  real dt_prof_cas(nlev_cas),ht_prof_cas(nlev_cas),vt_prof_cas(nlev_cas)
1087  real dtrad_prof_cas(nlev_cas)
1088  real dq_prof_cas(nlev_cas),hq_prof_cas(nlev_cas),vq_prof_cas(nlev_cas)
1089  real lat_prof_cas,sens_prof_cas,ts_prof_cas,ustar_prof_cas
1090  real uw_prof_cas(nlev_cas),vw_prof_cas(nlev_cas),q1_prof_cas(nlev_cas),q2_prof_cas(nlev_cas)
1091  ! local:
1092  integer it_cas1, it_cas2,k
1093  real timeit,time_cas1,time_cas2,frac
1094
1095
1096  print*,'Check time',day1,day_ju_ini_cas,day_deb+1,pdt_cas
1097
1098  ! On teste si la date du cas AMMA est correcte.
1099  ! C est pour memoire car en fait les fichiers .def
1100  ! sont censes etre corrects.
1101  ! A supprimer a terme (MPL 20150623)
1102  !     if ((forcing_type.eq.10).and.(1.eq.0)) then
1103  ! Check that initial day of the simulation consistent with AMMA case:
1104  !      if (annee_ref.ne.2006) then
1105  !       print*,'Pour AMMA, annee_ref doit etre 2006'
1106  !       print*,'Changer annee_ref dans run.def'
1107  !       stop
1108  !      endif
1109  !      if (annee_ref.eq.2006 .and. day1.lt.day_cas) then
1110  !       print*,'AMMA a debute le 10 juillet 2006',day1,day_cas
1111  !       print*,'Changer dayref dans run.def'
1112  !       stop
1113  !      endif
1114  !      if (annee_ref.eq.2006 .and. day1.gt.day_cas+1) then
1115  !       print*,'AMMA a fini le 11 juillet'
1116  !       print*,'Changer dayref ou nday dans run.def'
1117  !       stop
1118  !      endif
1119  !      endif
1120
1121  ! Determine timestep relative to the 1st day:
1122  !       timeit=(day-day1)*86400.
1123  !       if (annee_ref.eq.1992) then
1124  !        timeit=(day-day_cas)*86400.
1125  !       else
1126  !        timeit=(day+61.-1.)*86400. ! 61 days between Nov01 and Dec31 1992
1127  !       endif
1128  timeit=(day-day_ju_ini_cas)*86400
1129  !print *,'day=',day
1130  !print *,'day_ju_ini_cas=',day_ju_ini_cas
1131  !print *,'pdt_cas=',pdt_cas
1132  !print *,'timeit=',timeit
1133  !print *,'nt_cas=',nt_cas
1134
1135  ! Determine the closest observation times:
1136  !       it_cas1=INT(timeit/pdt_cas)+1
1137  !       it_cas2=it_cas1 + 1
1138  !       time_cas1=(it_cas1-1)*pdt_cas
1139  !       time_cas2=(it_cas2-1)*pdt_cas
1140
1141  it_cas1=INT(timeit/pdt_cas)+1
1142  IF (it_cas1 .EQ. nt_cas) THEN
1143     it_cas2=it_cas1
1144  ELSE
1145     it_cas2=it_cas1 + 1
1146  ENDIF
1147  time_cas1=(it_cas1-1)*pdt_cas
1148  time_cas2=(it_cas2-1)*pdt_cas
1149  !print *,'it_cas1,it_cas2,time_cas1,time_cas2=',it_cas1,it_cas2,time_cas1,time_cas2
1150
1151  if (it_cas1 .gt. nt_cas) then
1152     write(*,*) 'PB-stop: day, day_ju_ini_cas,it_cas1, it_cas2, timeit: '            &
1153          ,day,day_ju_ini_cas,it_cas1,it_cas2,timeit
1154     stop
1155  endif
1156
1157  ! time interpolation:
1158  IF (it_cas1 .EQ. it_cas2) THEN
1159     frac=0.
1160  ELSE
1161     frac=(time_cas2-timeit)/(time_cas2-time_cas1)
1162     frac=max(frac,0.0)
1163  ENDIF
1164
1165  lat_prof_cas = lat_cas(it_cas2)                                       &
1166       -frac*(lat_cas(it_cas2)-lat_cas(it_cas1))
1167  sens_prof_cas = sens_cas(it_cas2)                                     &
1168       -frac*(sens_cas(it_cas2)-sens_cas(it_cas1))
1169  ts_prof_cas = ts_cas(it_cas2)                                         &
1170       -frac*(ts_cas(it_cas2)-ts_cas(it_cas1))
1171  ustar_prof_cas = ustar_cas(it_cas2)                                   &
1172       -frac*(ustar_cas(it_cas2)-ustar_cas(it_cas1))
1173
1174  do k=1,nlev_cas
1175     plev_prof_cas(k) = plev_cas(k,it_cas2)                               &
1176          -frac*(plev_cas(k,it_cas2)-plev_cas(k,it_cas1))
1177     t_prof_cas(k) = t_cas(k,it_cas2)                               &
1178          -frac*(t_cas(k,it_cas2)-t_cas(k,it_cas1))
1179     q_prof_cas(k) = q_cas(k,it_cas2)                               &
1180          -frac*(q_cas(k,it_cas2)-q_cas(k,it_cas1))
1181     u_prof_cas(k) = u_cas(k,it_cas2)                               &
1182          -frac*(u_cas(k,it_cas2)-u_cas(k,it_cas1))
1183     v_prof_cas(k) = v_cas(k,it_cas2)                               &
1184          -frac*(v_cas(k,it_cas2)-v_cas(k,it_cas1))
1185     ug_prof_cas(k) = ug_cas(k,it_cas2)                               &
1186          -frac*(ug_cas(k,it_cas2)-ug_cas(k,it_cas1))
1187     vg_prof_cas(k) = vg_cas(k,it_cas2)                               &
1188          -frac*(vg_cas(k,it_cas2)-vg_cas(k,it_cas1))
1189     vitw_prof_cas(k) = vitw_cas(k,it_cas2)                               &
1190          -frac*(vitw_cas(k,it_cas2)-vitw_cas(k,it_cas1))
1191     du_prof_cas(k) = du_cas(k,it_cas2)                                   &
1192          -frac*(du_cas(k,it_cas2)-du_cas(k,it_cas1))
1193     hu_prof_cas(k) = hu_cas(k,it_cas2)                                   &
1194          -frac*(hu_cas(k,it_cas2)-hu_cas(k,it_cas1))
1195     vu_prof_cas(k) = vu_cas(k,it_cas2)                                   &
1196          -frac*(vu_cas(k,it_cas2)-vu_cas(k,it_cas1))
1197     dv_prof_cas(k) = dv_cas(k,it_cas2)                                   &
1198          -frac*(dv_cas(k,it_cas2)-dv_cas(k,it_cas1))
1199     hv_prof_cas(k) = hv_cas(k,it_cas2)                                   &
1200          -frac*(hv_cas(k,it_cas2)-hv_cas(k,it_cas1))
1201     vv_prof_cas(k) = vv_cas(k,it_cas2)                                   &
1202          -frac*(vv_cas(k,it_cas2)-vv_cas(k,it_cas1))
1203     dt_prof_cas(k) = dt_cas(k,it_cas2)                                   &
1204          -frac*(dt_cas(k,it_cas2)-dt_cas(k,it_cas1))
1205     ht_prof_cas(k) = ht_cas(k,it_cas2)                                   &
1206          -frac*(ht_cas(k,it_cas2)-ht_cas(k,it_cas1))
1207     vt_prof_cas(k) = vt_cas(k,it_cas2)                                   &
1208          -frac*(vt_cas(k,it_cas2)-vt_cas(k,it_cas1))
1209     dtrad_prof_cas(k) = dtrad_cas(k,it_cas2)                                   &
1210          -frac*(dtrad_cas(k,it_cas2)-dtrad_cas(k,it_cas1))
1211     dq_prof_cas(k) = dq_cas(k,it_cas2)                                   &
1212          -frac*(dq_cas(k,it_cas2)-dq_cas(k,it_cas1))
1213     hq_prof_cas(k) = hq_cas(k,it_cas2)                                   &
1214          -frac*(hq_cas(k,it_cas2)-hq_cas(k,it_cas1))
1215     vq_prof_cas(k) = vq_cas(k,it_cas2)                                   &
1216          -frac*(vq_cas(k,it_cas2)-vq_cas(k,it_cas1))
1217     uw_prof_cas(k) = uw_cas(k,it_cas2)                                   &
1218          -frac*(uw_cas(k,it_cas2)-uw_cas(k,it_cas1))
1219     vw_prof_cas(k) = vw_cas(k,it_cas2)                                   &
1220          -frac*(vw_cas(k,it_cas2)-vw_cas(k,it_cas1))
1221     q1_prof_cas(k) = q1_cas(k,it_cas2)                                   &
1222          -frac*(q1_cas(k,it_cas2)-q1_cas(k,it_cas1))
1223     q2_prof_cas(k) = q2_cas(k,it_cas2)                                   &
1224          -frac*(q2_cas(k,it_cas2)-q2_cas(k,it_cas1))
1225  enddo
1226
1227  return
1228END SUBROUTINE interp_case_time2
1229
1230!**********************************************************************************************
1231SUBROUTINE interp2_case_time(day,day1,annee_ref                           &
1232     !    &         ,year_cas,day_cas,nt_cas,pdt_forc,nlev_cas                         &
1233     ,nt_cas,nlev_cas                                                   &
1234     ,ts_cas,ps_cas,plev_cas,t_cas,theta_cas,thv_cas,thl_cas            &
1235     ,qv_cas,ql_cas,qi_cas,u_cas,v_cas                                  &
1236     ,ug_cas,vg_cas,vitw_cas,omega_cas,du_cas,hu_cas,vu_cas             &
1237     ,dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dtrad_cas               &
1238     ,dq_cas,hq_cas,vq_cas,dth_cas,hth_cas,vth_cas                      &
1239     ,lat_cas,sens_cas,ustar_cas                                        &
1240     ,uw_cas,vw_cas,q1_cas,q2_cas,tke_cas                               &
1241     !
1242     ,ts_prof_cas,plev_prof_cas,t_prof_cas,theta_prof_cas               &
1243     ,thv_prof_cas,thl_prof_cas,qv_prof_cas,ql_prof_cas,qi_prof_cas     &
1244     ,u_prof_cas,v_prof_cas,ug_prof_cas,vg_prof_cas                     &
1245     ,vitw_prof_cas,omega_prof_cas,du_prof_cas,hu_prof_cas,vu_prof_cas  &
1246     ,dv_prof_cas,hv_prof_cas,vv_prof_cas,dt_prof_cas                   &
1247     ,ht_prof_cas,vt_prof_cas,dtrad_prof_cas,dq_prof_cas                &
1248     ,hq_prof_cas,vq_prof_cas,dth_prof_cas,hth_prof_cas,vth_prof_cas    &
1249     ,lat_prof_cas,sens_prof_cas                                        &
1250     ,ustar_prof_cas,uw_prof_cas,vw_prof_cas,q1_prof_cas,q2_prof_cas,tke_prof_cas)
1251
1252  USE compar1d_mod_h
1253  USE date_cas_mod_h
1254  implicit none
1255
1256  !---------------------------------------------------------------------------------------
1257  ! Time interpolation of a 2D field to the timestep corresponding to day
1258  !
1259  ! day: current julian day (e.g. 717538.2)
1260  ! day1: first day of the simulation
1261  ! nt_cas: total nb of data in the forcing
1262  ! pdt_cas: total time interval (in sec) between 2 forcing data
1263  !---------------------------------------------------------------------------------------
1264
1265  ! inputs:
1266  integer annee_ref
1267  integer nt_cas,nlev_cas
1268  real day, day1,day_cas
1269  real ts_cas(nt_cas),ps_cas(nt_cas)
1270  real plev_cas(nlev_cas,nt_cas)
1271  real t_cas(nlev_cas,nt_cas),theta_cas(nlev_cas,nt_cas),thv_cas(nlev_cas,nt_cas),thl_cas(nlev_cas,nt_cas)
1272  real qv_cas(nlev_cas,nt_cas),ql_cas(nlev_cas,nt_cas),qi_cas(nlev_cas,nt_cas)
1273  real u_cas(nlev_cas,nt_cas),v_cas(nlev_cas,nt_cas)
1274  real ug_cas(nlev_cas,nt_cas),vg_cas(nlev_cas,nt_cas)
1275  real vitw_cas(nlev_cas,nt_cas),omega_cas(nlev_cas,nt_cas)
1276  real du_cas(nlev_cas,nt_cas),hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas)
1277  real dv_cas(nlev_cas,nt_cas),hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas)
1278  real dt_cas(nlev_cas,nt_cas),ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas)
1279  real dth_cas(nlev_cas,nt_cas),hth_cas(nlev_cas,nt_cas),vth_cas(nlev_cas,nt_cas)
1280  real dtrad_cas(nlev_cas,nt_cas)
1281  real dq_cas(nlev_cas,nt_cas),hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas)
1282  real lat_cas(nt_cas),sens_cas(nt_cas),tke_cas(nt_cas)
1283  real ustar_cas(nt_cas),uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas)
1284  real q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas)
1285
1286  ! outputs:
1287  real plev_prof_cas(nlev_cas)
1288  real t_prof_cas(nlev_cas),theta_prof_cas(nlev_cas),thl_prof_cas(nlev_cas),thv_prof_cas(nlev_cas)
1289  real qv_prof_cas(nlev_cas),ql_prof_cas(nlev_cas),qi_prof_cas(nlev_cas)
1290  real u_prof_cas(nlev_cas),v_prof_cas(nlev_cas)
1291  real ug_prof_cas(nlev_cas),vg_prof_cas(nlev_cas)
1292  real vitw_prof_cas(nlev_cas),omega_prof_cas(nlev_cas)
1293  real du_prof_cas(nlev_cas),hu_prof_cas(nlev_cas),vu_prof_cas(nlev_cas)
1294  real dv_prof_cas(nlev_cas),hv_prof_cas(nlev_cas),vv_prof_cas(nlev_cas)
1295  real dt_prof_cas(nlev_cas),ht_prof_cas(nlev_cas),vt_prof_cas(nlev_cas)
1296  real dth_prof_cas(nlev_cas),hth_prof_cas(nlev_cas),vth_prof_cas(nlev_cas)
1297  real dtrad_prof_cas(nlev_cas)
1298  real dq_prof_cas(nlev_cas),hq_prof_cas(nlev_cas),vq_prof_cas(nlev_cas)
1299  real lat_prof_cas,sens_prof_cas,tke_prof_cas,ts_prof_cas,ustar_prof_cas
1300  real uw_prof_cas(nlev_cas),vw_prof_cas(nlev_cas),q1_prof_cas(nlev_cas),q2_prof_cas(nlev_cas)
1301  ! local:
1302  integer it_cas1, it_cas2,k
1303  real timeit,time_cas1,time_cas2,frac
1304
1305
1306  print*,'Check time',day1,day_ju_ini_cas,day_deb+1,pdt_cas
1307  !       do k=1,nlev_cas
1308  !       print*,'debut de interp2_case_time, plev_cas=',k,plev_cas(k,1)
1309  !       enddo
1310
1311  ! On teste si la date du cas AMMA est correcte.
1312  ! C est pour memoire car en fait les fichiers .def
1313  ! sont censes etre corrects.
1314  ! A supprimer a terme (MPL 20150623)
1315  !     if ((forcing_type.eq.10).and.(1.eq.0)) then
1316  ! Check that initial day of the simulation consistent with AMMA case:
1317  !      if (annee_ref.ne.2006) then
1318  !       print*,'Pour AMMA, annee_ref doit etre 2006'
1319  !       print*,'Changer annee_ref dans run.def'
1320  !       stop
1321  !      endif
1322  !      if (annee_ref.eq.2006 .and. day1.lt.day_cas) then
1323  !       print*,'AMMA a debute le 10 juillet 2006',day1,day_cas
1324  !       print*,'Changer dayref dans run.def'
1325  !       stop
1326  !      endif
1327  !      if (annee_ref.eq.2006 .and. day1.gt.day_cas+1) then
1328  !       print*,'AMMA a fini le 11 juillet'
1329  !       print*,'Changer dayref ou nday dans run.def'
1330  !       stop
1331  !      endif
1332  !      endif
1333
1334  ! Determine timestep relative to the 1st day:
1335  !       timeit=(day-day1)*86400.
1336  !       if (annee_ref.eq.1992) then
1337  !        timeit=(day-day_cas)*86400.
1338  !       else
1339  !        timeit=(day+61.-1.)*86400. ! 61 days between Nov01 and Dec31 1992
1340  !       endif
1341  timeit=(day-day_ju_ini_cas)*86400
1342  !print *,'day=',day
1343  !print *,'day_ju_ini_cas=',day_ju_ini_cas
1344  !print *,'pdt_cas=',pdt_cas
1345  !print *,'timeit=',timeit
1346  !print *,'nt_cas=',nt_cas
1347
1348  ! Determine the closest observation times:
1349  !       it_cas1=INT(timeit/pdt_cas)+1
1350  !       it_cas2=it_cas1 + 1
1351  !       time_cas1=(it_cas1-1)*pdt_cas
1352  !       time_cas2=(it_cas2-1)*pdt_cas
1353
1354  it_cas1=INT(timeit/pdt_cas)+1
1355  IF (it_cas1 .EQ. nt_cas) THEN
1356     it_cas2=it_cas1
1357  ELSE
1358     it_cas2=it_cas1 + 1
1359  ENDIF
1360  time_cas1=(it_cas1-1)*pdt_cas
1361  time_cas2=(it_cas2-1)*pdt_cas
1362  !print *,'timeit,pdt_cas,nt_cas=',timeit,pdt_cas,nt_cas
1363  !print *,'it_cas1,it_cas2,time_cas1,time_cas2=',it_cas1,it_cas2,time_cas1,time_cas2
1364
1365  if (it_cas1 .gt. nt_cas) then
1366     write(*,*) 'PB-stop: day, day_ju_ini_cas,it_cas1, it_cas2, timeit: '            &
1367          ,day,day_ju_ini_cas,it_cas1,it_cas2,timeit
1368     stop
1369  endif
1370
1371  ! time interpolation:
1372  IF (it_cas1 .EQ. it_cas2) THEN
1373     frac=0.
1374  ELSE
1375     frac=(time_cas2-timeit)/(time_cas2-time_cas1)
1376     frac=max(frac,0.0)
1377  ENDIF
1378
1379  lat_prof_cas = lat_cas(it_cas2)                                   &
1380       -frac*(lat_cas(it_cas2)-lat_cas(it_cas1))
1381  sens_prof_cas = sens_cas(it_cas2)                                 &
1382       -frac*(sens_cas(it_cas2)-sens_cas(it_cas1))
1383  tke_prof_cas = tke_cas(it_cas2)                                   &
1384       -frac*(tke_cas(it_cas2)-tke_cas(it_cas1))
1385  ts_prof_cas = ts_cas(it_cas2)                                     &
1386       -frac*(ts_cas(it_cas2)-ts_cas(it_cas1))
1387  ustar_prof_cas = ustar_cas(it_cas2)                               &
1388       -frac*(ustar_cas(it_cas2)-ustar_cas(it_cas1))
1389
1390  do k=1,nlev_cas
1391     plev_prof_cas(k) = plev_cas(k,it_cas2)                           &     
1392          -frac*(plev_cas(k,it_cas2)-plev_cas(k,it_cas1))
1393     t_prof_cas(k) = t_cas(k,it_cas2)                                 &       
1394          -frac*(t_cas(k,it_cas2)-t_cas(k,it_cas1))
1395     !print *,'k,frac,plev_cas1,plev_cas2=',k,frac,plev_cas(k,it_cas1),plev_cas(k,it_cas2)
1396     theta_prof_cas(k) = theta_cas(k,it_cas2)                         &                     
1397          -frac*(theta_cas(k,it_cas2)-theta_cas(k,it_cas1))
1398     thv_prof_cas(k) = thv_cas(k,it_cas2)                             &         
1399          -frac*(thv_cas(k,it_cas2)-thv_cas(k,it_cas1))
1400     thl_prof_cas(k) = thl_cas(k,it_cas2)                             &             
1401          -frac*(thl_cas(k,it_cas2)-thl_cas(k,it_cas1))
1402     qv_prof_cas(k) = qv_cas(k,it_cas2)                               &
1403          -frac*(qv_cas(k,it_cas2)-qv_cas(k,it_cas1))
1404     ql_prof_cas(k) = ql_cas(k,it_cas2)                               &
1405          -frac*(ql_cas(k,it_cas2)-ql_cas(k,it_cas1))
1406     qi_prof_cas(k) = qi_cas(k,it_cas2)                               &
1407          -frac*(qi_cas(k,it_cas2)-qi_cas(k,it_cas1))
1408     u_prof_cas(k) = u_cas(k,it_cas2)                                 &
1409          -frac*(u_cas(k,it_cas2)-u_cas(k,it_cas1))
1410     v_prof_cas(k) = v_cas(k,it_cas2)                                 &
1411          -frac*(v_cas(k,it_cas2)-v_cas(k,it_cas1))
1412     ug_prof_cas(k) = ug_cas(k,it_cas2)                               &
1413          -frac*(ug_cas(k,it_cas2)-ug_cas(k,it_cas1))
1414     vg_prof_cas(k) = vg_cas(k,it_cas2)                               &
1415          -frac*(vg_cas(k,it_cas2)-vg_cas(k,it_cas1))
1416     vitw_prof_cas(k) = vitw_cas(k,it_cas2)                           &
1417          -frac*(vitw_cas(k,it_cas2)-vitw_cas(k,it_cas1))
1418     omega_prof_cas(k) = omega_cas(k,it_cas2)                         &
1419          -frac*(omega_cas(k,it_cas2)-omega_cas(k,it_cas1))
1420     du_prof_cas(k) = du_cas(k,it_cas2)                               &
1421          -frac*(du_cas(k,it_cas2)-du_cas(k,it_cas1))
1422     hu_prof_cas(k) = hu_cas(k,it_cas2)                               &
1423          -frac*(hu_cas(k,it_cas2)-hu_cas(k,it_cas1))
1424     vu_prof_cas(k) = vu_cas(k,it_cas2)                               &
1425          -frac*(vu_cas(k,it_cas2)-vu_cas(k,it_cas1))
1426     dv_prof_cas(k) = dv_cas(k,it_cas2)                               &
1427          -frac*(dv_cas(k,it_cas2)-dv_cas(k,it_cas1))
1428     hv_prof_cas(k) = hv_cas(k,it_cas2)                               &
1429          -frac*(hv_cas(k,it_cas2)-hv_cas(k,it_cas1))
1430     vv_prof_cas(k) = vv_cas(k,it_cas2)                               &
1431          -frac*(vv_cas(k,it_cas2)-vv_cas(k,it_cas1))
1432     dt_prof_cas(k) = dt_cas(k,it_cas2)                               &
1433          -frac*(dt_cas(k,it_cas2)-dt_cas(k,it_cas1))
1434     ht_prof_cas(k) = ht_cas(k,it_cas2)                               &
1435          -frac*(ht_cas(k,it_cas2)-ht_cas(k,it_cas1))
1436     vt_prof_cas(k) = vt_cas(k,it_cas2)                               &
1437          -frac*(vt_cas(k,it_cas2)-vt_cas(k,it_cas1))
1438     dth_prof_cas(k) = dth_cas(k,it_cas2)                             &
1439          -frac*(dth_cas(k,it_cas2)-dth_cas(k,it_cas1))
1440     hth_prof_cas(k) = hth_cas(k,it_cas2)                             &
1441          -frac*(hth_cas(k,it_cas2)-hth_cas(k,it_cas1))
1442     vth_prof_cas(k) = vth_cas(k,it_cas2)                             &
1443          -frac*(vth_cas(k,it_cas2)-vth_cas(k,it_cas1))
1444     dtrad_prof_cas(k) = dtrad_cas(k,it_cas2)                         &
1445          -frac*(dtrad_cas(k,it_cas2)-dtrad_cas(k,it_cas1))
1446     dq_prof_cas(k) = dq_cas(k,it_cas2)                               &
1447          -frac*(dq_cas(k,it_cas2)-dq_cas(k,it_cas1))
1448     hq_prof_cas(k) = hq_cas(k,it_cas2)                               &
1449          -frac*(hq_cas(k,it_cas2)-hq_cas(k,it_cas1))
1450     vq_prof_cas(k) = vq_cas(k,it_cas2)                               &
1451          -frac*(vq_cas(k,it_cas2)-vq_cas(k,it_cas1))
1452     uw_prof_cas(k) = uw_cas(k,it_cas2)                                &
1453          -frac*(uw_cas(k,it_cas2)-uw_cas(k,it_cas1))
1454     vw_prof_cas(k) = vw_cas(k,it_cas2)                                &
1455          -frac*(vw_cas(k,it_cas2)-vw_cas(k,it_cas1))
1456     q1_prof_cas(k) = q1_cas(k,it_cas2)                                &
1457          -frac*(q1_cas(k,it_cas2)-q1_cas(k,it_cas1))
1458     q2_prof_cas(k) = q2_cas(k,it_cas2)                                &
1459          -frac*(q2_cas(k,it_cas2)-q2_cas(k,it_cas1))
1460  enddo
1461
1462  return
1463END SUBROUTINE interp2_case_time
1464
1465!**********************************************************************************************
1466
Note: See TracBrowser for help on using the repository browser.