source: LMDZ6/trunk/libf/phylmd/dyn1d/mod_1D_cases_read_std.F90 @ 3781

Last change on this file since 3781 was 3781, checked in by evignon, 4 years ago

Initialisation de la TKE pour les cas 1D (important pour GABLS1), Etienne

File size: 53.0 KB
Line 
1!
2! $Id: mod_1D_cases_read.F90 2373 2015-10-13 17:28:01Z jyg $
3!
4MODULE mod_1D_cases_read_std
5
6!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
7!Declarations specifiques au cas standard
8        character*80 :: fich_cas
9! Discr?tisation
10        integer nlev_cas, nt_cas
11
12
13!profils environnementaux
14        real, allocatable::  plev_cas(:,:),plevh_cas(:)
15        real, allocatable::  ap_cas(:),bp_cas(:)
16
17        real, allocatable::  z_cas(:,:),zh_cas(:)
18        real, allocatable::  t_cas(:,:),q_cas(:,:),qv_cas(:,:),ql_cas(:,:),qi_cas(:,:),rh_cas(:,:)
19        real, allocatable::  th_cas(:,:),thv_cas(:,:),thl_cas(:,:),rv_cas(:,:)
20        real, allocatable::  u_cas(:,:),v_cas(:,:),vitw_cas(:,:),omega_cas(:,:),tke_cas(:,:)
21
22!forcing
23        real, allocatable::  ht_cas(:,:),vt_cas(:,:),dt_cas(:,:),dtrad_cas(:,:)
24        real, allocatable::  hth_cas(:,:),vth_cas(:,:),dth_cas(:,:)
25        real, allocatable::  hq_cas(:,:),vq_cas(:,:),dq_cas(:,:)
26        real, allocatable::  hr_cas(:,:),vr_cas(:,:),dr_cas(:,:)
27        real, allocatable::  hu_cas(:,:),vu_cas(:,:),du_cas(:,:)
28        real, allocatable::  hv_cas(:,:),vv_cas(:,:),dv_cas(:,:)
29        real, allocatable::  ug_cas(:,:),vg_cas(:,:)
30        real, allocatable::  temp_nudg_cas(:,:),qv_nudg_cas(:,:),u_nudg_cas(:,:),v_nudg_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(:,:),tkes_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::  tke_prof_cas(:)
51        real, allocatable::  ug_prof_cas(:)
52        real, allocatable::  vg_prof_cas(:)
53        real, allocatable::  temp_nudg_prof_cas(:),qv_nudg_prof_cas(:),u_nudg_prof_cas(:),v_nudg_prof_cas(:)
54        real, allocatable::  ht_prof_cas(:)
55        real, allocatable::  hth_prof_cas(:)
56        real, allocatable::  hq_prof_cas(:)
57        real, allocatable::  vt_prof_cas(:)
58        real, allocatable::  vth_prof_cas(:)
59        real, allocatable::  vq_prof_cas(:)
60        real, allocatable::  dt_prof_cas(:)
61        real, allocatable::  dth_prof_cas(:)
62        real, allocatable::  dtrad_prof_cas(:)
63        real, allocatable::  dq_prof_cas(:)
64        real, allocatable::  hu_prof_cas(:)
65        real, allocatable::  hv_prof_cas(:)
66        real, allocatable::  vu_prof_cas(:)
67        real, allocatable::  vv_prof_cas(:)
68        real, allocatable::  du_prof_cas(:)
69        real, allocatable::  dv_prof_cas(:)
70        real, allocatable::  uw_prof_cas(:)
71        real, allocatable::  vw_prof_cas(:)
72        real, allocatable::  q1_prof_cas(:)
73        real, allocatable::  q2_prof_cas(:)
74
75
76        real lat_prof_cas,sens_prof_cas,ts_prof_cas,ps_prof_cas,ustar_prof_cas,tkes_prof_cas
77        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
78     
79
80
81CONTAINS
82
83
84!**********************************************************************************************
85SUBROUTINE read_SCM_cas
86      implicit none
87
88#include "netcdf.inc"
89#include "date_cas.h"
90
91      INTEGER nid,rid,ierr
92      INTEGER ii,jj,timeid
93      REAL, ALLOCATABLE :: time_val(:)
94
95      print*,'ON EST VRAIMENT DASN MOD_1D_CASES_READ_STD'
96      fich_cas='cas.nc'
97      print*,'fich_cas ',fich_cas
98      ierr = NF_OPEN(fich_cas,NF_NOWRITE,nid)
99      print*,'fich_cas,NF_NOWRITE,nid ',fich_cas,NF_NOWRITE,nid
100      if (ierr.NE.NF_NOERR) then
101         write(*,*) 'ERROR: GROS Pb opening forcings nc file '
102         write(*,*) NF_STRERROR(ierr)
103         stop ""
104      endif
105!.......................................................................
106      ierr=NF_INQ_DIMID(nid,'lat',rid)
107      IF (ierr.NE.NF_NOERR) THEN
108         print*, 'Oh probleme lecture dimension lat'
109      ENDIF
110      ierr=NF_INQ_DIMLEN(nid,rid,ii)
111      print*,'OK1 read2: nid,rid,lat',nid,rid,ii
112!.......................................................................
113      ierr=NF_INQ_DIMID(nid,'lon',rid)
114      IF (ierr.NE.NF_NOERR) THEN
115         print*, 'Oh probleme lecture dimension lon'
116      ENDIF
117      ierr=NF_INQ_DIMLEN(nid,rid,jj)
118      print*,'OK2 read2: nid,rid,lat',nid,rid,jj
119!.......................................................................
120      ierr=NF_INQ_DIMID(nid,'lev',rid)
121      IF (ierr.NE.NF_NOERR) THEN
122         print*, 'Oh probleme lecture dimension nlev'
123      ENDIF
124      ierr=NF_INQ_DIMLEN(nid,rid,nlev_cas)
125      print*,'OK3 read2: nid,rid,nlev_cas',nid,rid,nlev_cas
126      IF ( .NOT. ( nlev_cas > 10 .AND. nlev_cas < 200000 )) THEN
127              print*,'Valeur de nlev_cas peu probable'
128              STOP
129      ENDIF
130!.......................................................................
131      ierr=NF_INQ_DIMID(nid,'time',rid)
132      nt_cas=0
133      IF (ierr.NE.NF_NOERR) THEN
134        stop 'Oh probleme lecture dimension time'
135      ENDIF
136      ierr=NF_INQ_DIMLEN(nid,rid,nt_cas)
137      print*,'OK4 read2: nid,rid,nt_cas',nid,rid,nt_cas
138! Lecture de l'axe des temps
139      print*,'LECTURE DU TEMPS'
140      ierr=NF_INQ_VARID(nid,'time',timeid)
141         if(ierr/=NF_NOERR) then
142           print *,'Variable time manquante dans cas.nc:'
143           ierr=NF_NOERR
144         else
145                 allocate(time_val(nt_cas))
146#ifdef NC_DOUBLE
147         ierr = NF_GET_VAR_DOUBLE(nid,timeid,time_val)
148#else
149           ierr = NF_GET_VAR_REAL(nid,timeid,time_val)
150#endif
151           if(ierr/=NF_NOERR) then
152              print *,'Pb a la lecture de time cas.nc: '
153           endif
154   endif
155   IF (nt_cas>1) THEN
156           pdt_cas=time_val(2)-time_val(1)
157   ELSE
158           pdt_cas=0.
159   ENDIF
160
161
162!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
163!profils moyens:
164        allocate(plev_cas(nlev_cas,nt_cas),plevh_cas(nlev_cas+1))       
165        allocate(z_cas(nlev_cas,nt_cas),zh_cas(nlev_cas+1))
166        allocate(ap_cas(nlev_cas+1),bp_cas(nt_cas+1))
167        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), &
168             qi_cas(nlev_cas,nt_cas),rh_cas(nlev_cas,nt_cas))
169        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))
170        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))
171        allocate(tke_cas(nlev_cas,nt_cas))
172!forcing
173        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))
174        allocate(hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas),dq_cas(nlev_cas,nt_cas))
175        allocate(hth_cas(nlev_cas,nt_cas),vth_cas(nlev_cas,nt_cas),dth_cas(nlev_cas,nt_cas))
176        allocate(hr_cas(nlev_cas,nt_cas),vr_cas(nlev_cas,nt_cas),dr_cas(nlev_cas,nt_cas))
177        allocate(hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas),du_cas(nlev_cas,nt_cas))
178        allocate(hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas),dv_cas(nlev_cas,nt_cas))
179        allocate(ug_cas(nlev_cas,nt_cas),vg_cas(nlev_cas,nt_cas))
180        allocate(temp_nudg_cas(nlev_cas,nt_cas),qv_nudg_cas(nlev_cas,nt_cas))
181        allocate(u_nudg_cas(nlev_cas,nt_cas),v_nudg_cas(nlev_cas,nt_cas))
182        allocate(lat_cas(nt_cas),sens_cas(nt_cas),ts_cas(nt_cas),ps_cas(nt_cas),ustar_cas(nt_cas),tkes_cas(nt_cas))
183        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))
184
185
186
187!champs interpoles
188        allocate(plev_prof_cas(nlev_cas))
189        allocate(t_prof_cas(nlev_cas))
190        allocate(theta_prof_cas(nlev_cas))
191        allocate(thl_prof_cas(nlev_cas))
192        allocate(thv_prof_cas(nlev_cas))
193        allocate(q_prof_cas(nlev_cas))
194        allocate(qv_prof_cas(nlev_cas))
195        allocate(ql_prof_cas(nlev_cas))
196        allocate(qi_prof_cas(nlev_cas))
197        allocate(rh_prof_cas(nlev_cas))
198        allocate(rv_prof_cas(nlev_cas))
199        allocate(u_prof_cas(nlev_cas))
200        allocate(v_prof_cas(nlev_cas))
201        allocate(vitw_prof_cas(nlev_cas))
202        allocate(omega_prof_cas(nlev_cas))
203        allocate(tke_prof_cas(nlev_cas))
204        allocate(ug_prof_cas(nlev_cas))
205        allocate(vg_prof_cas(nlev_cas))
206        allocate(temp_nudg_prof_cas(nlev_cas),qv_nudg_prof_cas(nlev_cas))
207        allocate(u_nudg_prof_cas(nlev_cas),v_nudg_prof_cas(nlev_cas))
208        allocate(ht_prof_cas(nlev_cas))
209        allocate(hth_prof_cas(nlev_cas))
210        allocate(hq_prof_cas(nlev_cas))
211        allocate(hu_prof_cas(nlev_cas))
212        allocate(hv_prof_cas(nlev_cas))
213        allocate(vt_prof_cas(nlev_cas))
214        allocate(vth_prof_cas(nlev_cas))
215        allocate(vq_prof_cas(nlev_cas))
216        allocate(vu_prof_cas(nlev_cas))
217        allocate(vv_prof_cas(nlev_cas))
218        allocate(dt_prof_cas(nlev_cas))
219        allocate(dth_prof_cas(nlev_cas))
220        allocate(dtrad_prof_cas(nlev_cas))
221        allocate(dq_prof_cas(nlev_cas))
222        allocate(du_prof_cas(nlev_cas))
223        allocate(dv_prof_cas(nlev_cas))
224        allocate(uw_prof_cas(nlev_cas))
225        allocate(vw_prof_cas(nlev_cas))
226        allocate(q1_prof_cas(nlev_cas))
227        allocate(q2_prof_cas(nlev_cas))
228
229        print*,'Allocations OK'
230        CALL read_SCM (nid,nlev_cas,nt_cas,                                                                     &
231     &     ap_cas,bp_cas,z_cas,plev_cas,zh_cas,plevh_cas,t_cas,th_cas,thv_cas,thl_cas,qv_cas,                   &
232     &     ql_cas,qi_cas,rh_cas,rv_cas,u_cas,v_cas,vitw_cas,omega_cas,tke_cas,ug_cas,vg_cas,                            &
233     &     temp_nudg_cas,qv_nudg_cas,u_nudg_cas,v_nudg_cas,                                                     &
234     &     du_cas,hu_cas,vu_cas,                                                                                &
235     &     dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dq_cas,hq_cas,vq_cas,dth_cas,hth_cas,vth_cas,               &
236     &     dr_cas,hr_cas,vr_cas,dtrad_cas,sens_cas,lat_cas,ts_cas,ps_cas,ustar_cas,tkes_cas,                      &
237     &     uw_cas,vw_cas,q1_cas,q2_cas,orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,heat_rough, &
238     &     o3_cas,rugos_cas,clay_cas,sand_cas)
239        print*,'read_SCM cas OK'
240        do ii=1,nlev_cas
241        print*,'apres read2_SCM, plev_cas=',ii,plev_cas(ii,1)
242        !print*,'apres read_SCM, plev_cas=',ii,omega_cas(ii,nt_cas/2+1)
243        enddo
244
245
246END SUBROUTINE read_SCM_cas
247
248
249!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
250SUBROUTINE deallocate2_1D_cases
251!profils environnementaux:
252        deallocate(plev_cas,plevh_cas)
253       
254        deallocate(z_cas,zh_cas)
255        deallocate(ap_cas,bp_cas)
256        deallocate(t_cas,q_cas,qv_cas,ql_cas,qi_cas,rh_cas)
257        deallocate(th_cas,thl_cas,thv_cas,rv_cas)
258        deallocate(u_cas,v_cas,vitw_cas,omega_cas,tke_cas)
259       
260!forcing
261        deallocate(ht_cas,vt_cas,dt_cas,dtrad_cas)
262        deallocate(hq_cas,vq_cas,dq_cas)
263        deallocate(hth_cas,vth_cas,dth_cas)
264        deallocate(hr_cas,vr_cas,dr_cas)
265        deallocate(hu_cas,vu_cas,du_cas)
266        deallocate(hv_cas,vv_cas,dv_cas)
267        deallocate(ug_cas)
268        deallocate(vg_cas)
269        deallocate(lat_cas,sens_cas,ts_cas,ps_cas,ustar_cas,tkes_cas,uw_cas,vw_cas,q1_cas,q2_cas)
270
271!champs interpoles
272        deallocate(plev_prof_cas)
273        deallocate(t_prof_cas)
274        deallocate(theta_prof_cas)
275        deallocate(thl_prof_cas)
276        deallocate(thv_prof_cas)
277        deallocate(q_prof_cas)
278        deallocate(qv_prof_cas)
279        deallocate(ql_prof_cas)
280        deallocate(qi_prof_cas)
281        deallocate(rh_prof_cas)
282        deallocate(rv_prof_cas)
283        deallocate(u_prof_cas)
284        deallocate(v_prof_cas)
285        deallocate(vitw_prof_cas)
286        deallocate(omega_prof_cas)
287        deallocate(tke_prof_cas)
288        deallocate(ug_prof_cas)
289        deallocate(vg_prof_cas)
290        deallocate(temp_nudg_prof_cas,qv_nudg_prof_cas,u_nudg_prof_cas,v_nudg_prof_cas)
291        deallocate(ht_prof_cas)
292        deallocate(hq_prof_cas)
293        deallocate(hu_prof_cas)
294        deallocate(hv_prof_cas)
295        deallocate(vt_prof_cas)
296        deallocate(vq_prof_cas)
297        deallocate(vu_prof_cas)
298        deallocate(vv_prof_cas)
299        deallocate(dt_prof_cas)
300        deallocate(dtrad_prof_cas)
301        deallocate(dq_prof_cas)
302        deallocate(du_prof_cas)
303        deallocate(dv_prof_cas)
304        deallocate(t_prof_cas)
305        deallocate(u_prof_cas)
306        deallocate(v_prof_cas)
307        deallocate(uw_prof_cas)
308        deallocate(vw_prof_cas)
309        deallocate(q1_prof_cas)
310        deallocate(q2_prof_cas)
311
312END SUBROUTINE deallocate2_1D_cases
313
314
315!=====================================================================
316      SUBROUTINE read_SCM(nid,nlevel,ntime,                                       &
317     &     ap,bp,zz,pp,zzh,pph,temp,theta,thv,thl,qv,ql,qi,rh,rv,u,v,vitw,omega,tke,ug,vg,&
318     &     temp_nudg,qv_nudg,u_nudg,v_nudg,                                        &
319     &     du,hu,vu,dv,hv,vv,dt,ht,vt,dq,hq,vq,                                    &
320     &     dth,hth,vth,dr,hr,vr,dtrad,sens,flat,ts,ps,ustar,tkes,uw,vw,q1,q2,       &
321     &     orog_cas,albedo_cas,emiss_cas,t_skin_cas,q_skin_cas,mom_rough,          &
322     &     heat_rough,o3_cas,rugos_cas,clay_cas,sand_cas)
323
324!program reading forcing of the case study
325      implicit none
326#include "netcdf.inc"
327#include "compar1d.h"
328
329      integer ntime,nlevel,k,t
330
331      real ap(nlevel+1),bp(nlevel+1)
332      real zz(nlevel,ntime),zzh(nlevel+1)
333      real pp(nlevel,ntime),pph(nlevel+1)
334!profils initiaux
335      real temp0(nlevel),qv0(nlevel),ql0(nlevel),qi0(nlevel),u0(nlevel),v0(nlevel),tke0(nlevel)
336      real pp0(nlevel)   
337      real temp(nlevel,ntime),qv(nlevel,ntime),ql(nlevel,ntime),qi(nlevel,ntime),rh(nlevel,ntime)
338      real theta(nlevel,ntime),thv(nlevel,ntime),thl(nlevel,ntime),rv(nlevel,ntime)
339      real u(nlevel,ntime),v(nlevel,ntime),tkes(ntime)
340      real temp_nudg(nlevel,ntime),qv_nudg(nlevel,ntime),u_nudg(nlevel,ntime),v_nudg(nlevel,ntime)
341      real ug(nlevel,ntime),vg(nlevel,ntime)
342      real vitw(nlevel,ntime),omega(nlevel,ntime),tke(nlevel,ntime)
343      real du(nlevel,ntime),hu(nlevel,ntime),vu(nlevel,ntime)
344      real dv(nlevel,ntime),hv(nlevel,ntime),vv(nlevel,ntime)
345      real dt(nlevel,ntime),ht(nlevel,ntime),vt(nlevel,ntime)
346      real dtrad(nlevel,ntime)
347      real dq(nlevel,ntime),hq(nlevel,ntime),vq(nlevel,ntime)
348      real dth(nlevel,ntime),hth(nlevel,ntime),vth(nlevel,ntime),hthl(nlevel,ntime)
349      real dr(nlevel,ntime),hr(nlevel,ntime),vr(nlevel,ntime)
350      real flat(ntime),sens(ntime),ustar(ntime)
351      real uw(nlevel,ntime),vw(nlevel,ntime),q1(nlevel,ntime),q2(nlevel,ntime)
352      real ts(ntime),ps(ntime)
353      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
354      real apbp(nlevel+1),resul(nlevel,ntime),resul1(nlevel),resul2(ntime),resul3
355
356
357      integer nid, ierr,ierr1,ierr2,rid,i
358      integer nbvar3d
359      parameter(nbvar3d=74)
360      integer var3didin(nbvar3d),missing_var(nbvar3d)
361      character*13 name_var(1:nbvar3d)
362
363
364      data name_var/ &
365     ! coordonnees pression (n+1 niveaux) #4
366     & 'coor_par_a','coor_par_b','height_h','pressure_h',& ! #1-#4
367     ! coordonnees pression (n niveaux) #8
368     &'temp','qv','ql','qi','u','v','tke','pressure',& ! #5-#12
369     ! coordonnees pression + temps #42
370     &'w','omega','ug','vg','uadv','uadvh','uadvv','vadv','vadvh','vadvv','temp_adv','tadvh','tadvv',& !  #13 - #25
371     &'qv_adv','qvadvh','qvadvv','thadv','thadvh','thadvv','thladvh',                             & ! #26 - #32
372     & 'radv','radvh','radvv','radcool','q1','q2','ustress','vstress',                           & ! #33 - #40
373     & 'rh','temp_nudging','qv_nudging','u_nudging','v_nudging',                                       & ! #41-45
374     &'height_f','pressure_forc','tempt','theta','thv','thl','qvt','qlt','qit','rv','ut','vt',   & ! #46-58
375     ! coordonnees temps #12
376     &'tkes','sfc_sens_flx','sfc_lat_flx','ts','ps','ustar',&
377     &'orog','albedo','emiss','t_skin','q_skin','mom_rough','heat_rough',&
378     ! scalaires #4
379     &'o3','rugos','clay','sand'/
380
381!-----------------------------------------------------------------------
382! Checking availability of variable #i in the cas.nc file
383!     missing_var=1 if the variable is missing
384!-----------------------------------------------------------------------
385
386       do i=1,nbvar3d
387         missing_var(i)=0.
388         ierr=NF_INQ_VARID(nid,name_var(i),var3didin(i))
389         if(ierr/=NF_NOERR) then
390           print *,'Variable manquante dans cas.nc:',i,name_var(i)
391           ierr=NF_NOERR
392           missing_var(i)=1
393         else
394
395!-----------------------------------------------------------------------
396! Activating keys depending on the presence of specific variables in cas.nc
397!-----------------------------------------------------------------------
398if ( 1 == 1 ) THEN
399! A MODIFIER: il faudrait dire nudging_temp mais faut le declarer dans compar1d.h etc...       
400!           if ( name_var(i) == 'temp_nudging' .and. nint(nudging_t)==0) stop 'Nudging inconsistency temp'
401            if ( name_var(i) == 'qv_nudging' .and. nint(nudging_qv)==0) stop 'Nudging inconsistency qv'
402            if ( name_var(i) == 'u_nudging' .and. nint(nudging_u)==0) stop 'Nudging inconsistency u'
403            if ( name_var(i) == 'v_nudging' .and. nint(nudging_u)==0) stop 'Nudging inconsistency v'
404    ELSE
405             print*,'GUIDAGE : CONSISTENCY CHECK DEACTIVATED FOR TESTS of SANDU/REF'
406    ENDIF
407
408!-----------------------------------------------------------------------
409! Reading variables 1D (N+1) vertical variables (nlevelp1,lat,lon)
410!-----------------------------------------------------------------------
411           if(i.LE.4) then
412#ifdef NC_DOUBLE
413           ierr = NF_GET_VAR_DOUBLE(nid,var3didin(i),apbp)
414#else
415           ierr = NF_GET_VAR_REAL(nid,var3didin(i),apbp)
416#endif
417           print *,'read2_cas(apbp), on a lu ',i,name_var(i)
418           if(ierr/=NF_NOERR) then
419              print *,'Pb a la lecture de cas.nc: ',name_var(i)
420              stop "getvarup"
421           endif
422
423!-----------------------------------------------------------------------
424!  Reading 1D (N) vertical varialbes    (nlevel,lat,lon)   
425!-----------------------------------------------------------------------
426           else if(i.gt.4.and.i.LE.12) then 
427#ifdef NC_DOUBLE
428           ierr = NF_GET_VAR_DOUBLE(nid,var3didin(i),resul1)
429#else
430           ierr = NF_GET_VAR_REAL(nid,var3didin(i),resul1)
431#endif
432           print *,'read2_cas(resul1), on a lu ',i,name_var(i)
433           if(ierr/=NF_NOERR) then
434              print *,'Pb a la lecture de cas.nc: ',name_var(i)
435              stop "getvarup"
436           endif
437         print*,'Lecture de la variable #i ',i,name_var(i),minval(resul1),maxval(resul1)
438
439!-----------------------------------------------------------------------
440!  Reading 2D tim-vertical variables  (time,nlevel,lat,lon)
441!  TBD : seems to be the same as above.
442!-----------------------------------------------------------------------
443           else if(i.gt.12.and.i.LE.57) then
444#ifdef NC_DOUBLE
445           ierr = NF_GET_VAR_DOUBLE(nid,var3didin(i),resul)
446#else
447           ierr = NF_GET_VAR_REAL(nid,var3didin(i),resul)
448#endif
449           print *,'read2_cas(resul), on a lu ',i,name_var(i)
450           if(ierr/=NF_NOERR) then
451              print *,'Pb a la lecture de cas.nc: ',name_var(i)
452              stop "getvarup"
453           endif
454         print*,'Lecture de la variable #i ',i,name_var(i),minval(resul),maxval(resul)
455
456!-----------------------------------------------------------------------
457!  Reading 1D time variables (time,lat,lon)
458!-----------------------------------------------------------------------
459           else if (i.gt.57.and.i.LE.63) then
460#ifdef NC_DOUBLE
461           ierr = NF_GET_VAR_DOUBLE(nid,var3didin(i),resul2)
462#else
463           ierr = NF_GET_VAR_REAL(nid,var3didin(i),resul2)
464#endif
465           print *,'read2_cas(resul2), on a lu ',i,name_var(i)
466           if(ierr/=NF_NOERR) then
467              print *,'Pb a la lecture de cas.nc: ',name_var(i)
468              stop "getvarup"
469           endif
470         print*,'Lecture de la variable #i  ',i,name_var(i),minval(resul2),maxval(resul2)
471
472!-----------------------------------------------------------------------
473! Reading scalar variables (lat,lon)
474!-----------------------------------------------------------------------
475           else
476#ifdef NC_DOUBLE
477           ierr = NF_GET_VAR_DOUBLE(nid,var3didin(i),resul3)
478#else
479           ierr = NF_GET_VAR_REAL(nid,var3didin(i),resul3)
480#endif
481           print *,'read2_cas(resul3), on a lu ',i,name_var(i)
482           if(ierr/=NF_NOERR) then
483              print *,'Pb a la lecture de cas.nc: ',name_var(i)
484              stop "getvarup"
485           endif
486         print*,'Lecture de la variable #i ',i,name_var(i),resul3
487           endif
488         endif
489
490!-----------------------------------------------------------------------
491! Attributing variables
492!-----------------------------------------------------------------------
493         select case(i)
494         !case(1) ; ap=apbp       ! donnees indexees en nlevel+1
495         ! case(2) ; bp=apbp
496           case(3) ; zzh=apbp
497           case(4) ; pph=apbp
498           case(5) ; temp0=resul1    ! donnees initiales
499           case(6) ; qv0=resul1
500           case(7) ; ql0=resul1
501           case(8) ; qi0=resul1
502           case(9) ; u0=resul1
503           case(10) ; v0=resul1
504           case(11) ; tke0=resul1
505           case(12) ; pp0=resul1
506           case(13) ; vitw=resul    ! donnees indexees en nlevel,time
507           case(14) ; omega=resul
508           case(15) ; ug=resul
509           case(16) ; vg=resul
510           case(17) ; du=resul
511           case(18) ; hu=resul
512           case(19) ; vu=resul
513           case(20) ; dv=resul
514           case(21) ; hv=resul
515           case(22) ; vv=resul
516           case(23) ; dt=resul
517           case(24) ; ht=resul
518           case(25) ; vt=resul
519           case(26) ; dq=resul
520           case(27) ; hq=resul
521           case(28) ; vq=resul
522           case(29) ; dth=resul
523           case(30) ; hth=resul
524           case(31) ; vth=resul
525           case(32) ; hthl=resul
526           case(33) ; dr=resul
527           case(34) ; hr=resul
528           case(35) ; vr=resul
529           case(36) ; dtrad=resul
530           case(37) ; q1=resul
531           case(38) ; q2=resul
532           case(39) ; uw=resul
533           case(40) ; vw=resul
534           case(41) ; rh=resul
535           case(42) ; temp_nudg=resul
536           case(43) ; qv_nudg=resul
537           case(44) ; u_nudg=resul
538           case(45) ; v_nudg=resul
539           case(46) ; zz=resul      ! donnees en time,nlevel pour profil initial
540           case(47) ; pp=resul
541           case(48) ; temp=resul
542           case(49) ; theta=resul
543           case(50) ; thv=resul
544           case(51) ; thl=resul
545           case(52) ; qv=resul
546           case(53) ; ql=resul
547           case(54) ; qi=resul
548           case(55) ; rv=resul
549           case(56) ; u=resul
550           case(57) ; v=resul
551           case(58) ; tkes=resul2   ! donnees indexees en time
552           case(59) ; sens=resul2
553           case(60) ; flat=resul2
554           case(61) ; ts=resul2
555           case(62) ; ps=resul2
556           case(63) ; ustar=resul2
557           case(64) ; orog_cas=resul3      ! constantes
558           case(65) ; albedo_cas=resul3
559           case(66) ; emiss_cas=resul3
560           case(67) ; t_skin_cas=resul3
561           case(68) ; q_skin_cas=resul3
562           case(69) ; mom_rough=resul3
563           case(70) ; heat_rough=resul3
564           case(71) ; o3_cas=resul3       
565           case(72) ; rugos_cas=resul3
566           case(73) ; clay_cas=resul3
567           case(74) ; sand_cas=resul3
568         end select
569         resul=0.
570         resul1=0.
571         resul2=0.
572         resul3=0.
573       enddo
574         print*,'Lecture de la variable APRES ,sens ',minval(sens),maxval(sens)
575         print*,'Lecture de la variable APRES ,flat ',minval(flat),maxval(flat)
576
577!CR:ATTENTION EN ATTENTE DE REGLER LA QUESTION DU PAS DE TEMPS INITIAL
578       do t=1,ntime
579          do k=1,nlevel
580             temp(k,t)=temp0(k)
581             qv(k,t)=qv0(k)
582             ql(k,t)=ql0(k)
583             qi(k,t)=qi0(k)
584             u(k,t)=u0(k)
585             v(k,t)=v0(k)
586             tke(k,t)=tke0(k)
587          enddo
588       enddo
589       !!!! TRAVAIL : EN FONCTION DES DECISIONS SUR LA SPECIFICATION DE W
590       !!!omega=-vitw*pres*rg/(rd*temp)
591!-----------------------------------------------------------------------
592
593         return
594         END SUBROUTINE read_SCM
595!======================================================================
596
597!======================================================================
598
599!**********************************************************************************************
600
601!**********************************************************************************************
602        SUBROUTINE interp_case_time_std(day,day1,annee_ref                           &
603!    &         ,year_cas,day_cas,nt_cas,pdt_forc,nlev_cas                         &
604     &         ,nt_cas,nlev_cas                                                   &
605     &         ,ts_cas,ps_cas,plev_cas,t_cas,theta_cas,thv_cas,thl_cas            &
606     &         ,qv_cas,ql_cas,qi_cas,u_cas,v_cas                                  &
607     &         ,ug_cas,vg_cas,temp_nudg_cas,qv_nudg_cas,u_nudg_cas,v_nudg_cas     &
608     &         ,vitw_cas,omega_cas,tke_cas,du_cas,hu_cas,vu_cas             &
609     &         ,dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dtrad_cas               &
610     &         ,dq_cas,hq_cas,vq_cas,dth_cas,hth_cas,vth_cas                      &
611     &         ,lat_cas,sens_cas,ustar_cas                                        &
612     &         ,uw_cas,vw_cas,q1_cas,q2_cas,tkes_cas                               &
613!
614     &         ,ts_prof_cas,ps_prof_cas,plev_prof_cas,t_prof_cas,theta_prof_cas               &
615     &         ,thv_prof_cas,thl_prof_cas,qv_prof_cas,ql_prof_cas,qi_prof_cas     &
616     &         ,u_prof_cas,v_prof_cas,ug_prof_cas,vg_prof_cas                     &
617     &         ,temp_nudg_prof_cas,qv_nudg_prof_cas,u_nudg_prof_cas,v_nudg_prof_cas     &
618     &         ,vitw_prof_cas,omega_prof_cas,tke_prof_cas,du_prof_cas,hu_prof_cas,vu_prof_cas  &
619     &         ,dv_prof_cas,hv_prof_cas,vv_prof_cas,dt_prof_cas                   &
620     &         ,ht_prof_cas,vt_prof_cas,dtrad_prof_cas,dq_prof_cas                &
621     &         ,hq_prof_cas,vq_prof_cas,dth_prof_cas,hth_prof_cas,vth_prof_cas    &
622     &         ,lat_prof_cas,sens_prof_cas                                        &
623     &         ,ustar_prof_cas,uw_prof_cas,vw_prof_cas,q1_prof_cas,q2_prof_cas,tkes_prof_cas)
624
625
626
627
628
629
630       implicit none
631
632!---------------------------------------------------------------------------------------
633! Time interpolation of a 2D field to the timestep corresponding to day
634!
635! day: current julian day (e.g. 717538.2)
636! day1: first day of the simulation
637! nt_cas: total nb of data in the forcing
638! pdt_cas: total time interval (in sec) between 2 forcing data
639!---------------------------------------------------------------------------------------
640
641#include "compar1d.h"
642#include "date_cas.h"
643
644! inputs:
645        integer annee_ref
646        integer nt_cas,nlev_cas
647        real day, day1,day_cas
648        real ts_cas(nt_cas),ps_cas(nt_cas)
649        real plev_cas(nlev_cas,nt_cas)
650        real t_cas(nlev_cas,nt_cas),theta_cas(nlev_cas,nt_cas)
651        real thv_cas(nlev_cas,nt_cas), thl_cas(nlev_cas,nt_cas)
652        real qv_cas(nlev_cas,nt_cas),ql_cas(nlev_cas,nt_cas),qi_cas(nlev_cas,nt_cas)
653        real u_cas(nlev_cas,nt_cas),v_cas(nlev_cas,nt_cas)
654        real ug_cas(nlev_cas,nt_cas),vg_cas(nlev_cas,nt_cas)
655        real temp_nudg_cas(nlev_cas,nt_cas),qv_nudg_cas(nlev_cas,nt_cas)
656        real u_nudg_cas(nlev_cas,nt_cas),v_nudg_cas(nlev_cas,nt_cas)
657
658        real vitw_cas(nlev_cas,nt_cas),omega_cas(nlev_cas,nt_cas),tke_cas(nlev_cas,nt_cas)
659        real du_cas(nlev_cas,nt_cas),hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas)
660        real dv_cas(nlev_cas,nt_cas),hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas)
661        real dt_cas(nlev_cas,nt_cas),ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas)
662        real dth_cas(nlev_cas,nt_cas),hth_cas(nlev_cas,nt_cas),vth_cas(nlev_cas,nt_cas)
663        real dtrad_cas(nlev_cas,nt_cas)
664        real dq_cas(nlev_cas,nt_cas),hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas)
665        real lat_cas(nt_cas),sens_cas(nt_cas),tkes_cas(nt_cas)
666        real ustar_cas(nt_cas),uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas)
667        real q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas)
668
669! outputs:
670        real plev_prof_cas(nlev_cas)
671        real t_prof_cas(nlev_cas),theta_prof_cas(nlev_cas),thl_prof_cas(nlev_cas),thv_prof_cas(nlev_cas)
672        real qv_prof_cas(nlev_cas),ql_prof_cas(nlev_cas),qi_prof_cas(nlev_cas)
673        real u_prof_cas(nlev_cas),v_prof_cas(nlev_cas)
674        real ug_prof_cas(nlev_cas),vg_prof_cas(nlev_cas)
675        real temp_nudg_prof_cas(nlev_cas),qv_nudg_prof_cas(nlev_cas)
676        real u_nudg_prof_cas(nlev_cas),v_nudg_prof_cas(nlev_cas)
677
678        real vitw_prof_cas(nlev_cas),omega_prof_cas(nlev_cas),tke_prof_cas(nlev_cas)
679        real du_prof_cas(nlev_cas),hu_prof_cas(nlev_cas),vu_prof_cas(nlev_cas)
680        real dv_prof_cas(nlev_cas),hv_prof_cas(nlev_cas),vv_prof_cas(nlev_cas)
681        real dt_prof_cas(nlev_cas),ht_prof_cas(nlev_cas),vt_prof_cas(nlev_cas)
682        real dth_prof_cas(nlev_cas),hth_prof_cas(nlev_cas),vth_prof_cas(nlev_cas)
683        real dtrad_prof_cas(nlev_cas)
684        real dq_prof_cas(nlev_cas),hq_prof_cas(nlev_cas),vq_prof_cas(nlev_cas)
685        real lat_prof_cas,sens_prof_cas,tkes_prof_cas,ts_prof_cas,ps_prof_cas,ustar_prof_cas
686        real uw_prof_cas(nlev_cas),vw_prof_cas(nlev_cas),q1_prof_cas(nlev_cas),q2_prof_cas(nlev_cas)
687! local:
688        integer it_cas1, it_cas2,k
689        real timeit,time_cas1,time_cas2,frac
690
691
692        print*,'Check time',day1,day_ju_ini_cas,day_deb+1,pdt_cas
693!       do k=1,nlev_cas
694!       print*,'debut de interp2_case_time, plev_cas=',k,plev_cas(k,1)
695!       enddo
696
697! On teste si la date du cas AMMA est correcte.
698! C est pour memoire car en fait les fichiers .def
699! sont censes etre corrects.
700! A supprimer a terme (MPL 20150623)
701!     if ((forcing_type.eq.10).and.(1.eq.0)) then
702! Check that initial day of the simulation consistent with AMMA case:
703!      if (annee_ref.ne.2006) then
704!       print*,'Pour AMMA, annee_ref doit etre 2006'
705!       print*,'Changer annee_ref dans run.def'
706!       stop
707!      endif
708!      if (annee_ref.eq.2006 .and. day1.lt.day_cas) then
709!       print*,'AMMA a debute le 10 juillet 2006',day1,day_cas
710!       print*,'Changer dayref dans run.def'
711!       stop
712!      endif
713!      if (annee_ref.eq.2006 .and. day1.gt.day_cas+1) then
714!       print*,'AMMA a fini le 11 juillet'
715!       print*,'Changer dayref ou nday dans run.def'
716!       stop
717!      endif
718!      endif
719
720! Determine timestep relative to the 1st day:
721!       timeit=(day-day1)*86400.
722!       if (annee_ref.eq.1992) then
723!        timeit=(day-day_cas)*86400.
724!       else
725!        timeit=(day+61.-1.)*86400. ! 61 days between Nov01 and Dec31 1992
726!       endif
727      timeit=(day-day_ju_ini_cas)*86400
728      print *,'day=',day
729      print *,'day_ju_ini_cas=',day_ju_ini_cas
730      print *,'pdt_cas=',pdt_cas
731      print *,'timeit=',timeit
732      print *,'nt_cas=',nt_cas
733
734! Determine the closest observation times:
735!       it_cas1=INT(timeit/pdt_cas)+1
736!       it_cas2=it_cas1 + 1
737!       time_cas1=(it_cas1-1)*pdt_cas
738!       time_cas2=(it_cas2-1)*pdt_cas
739
740       it_cas1=INT(timeit/pdt_cas)+1
741       IF (it_cas1 .EQ. nt_cas) THEN
742       it_cas2=it_cas1
743       ELSE
744       it_cas2=it_cas1 + 1
745       ENDIF
746       time_cas1=(it_cas1-1)*pdt_cas
747       time_cas2=(it_cas2-1)*pdt_cas
748!     print *,'timeit,pdt_cas,nt_cas=',timeit,pdt_cas,nt_cas
749!     print *,'it_cas1,it_cas2,time_cas1,time_cas2=',it_cas1,it_cas2,time_cas1,time_cas2
750
751       if (it_cas1 .gt. nt_cas) then
752        write(*,*) 'PB-stop: day, day_ju_ini_cas,it_cas1, it_cas2, timeit: '            &
753     &        ,day,day_ju_ini_cas,it_cas1,it_cas2,timeit
754        stop
755       endif
756
757! time interpolation:
758       IF (it_cas1 .EQ. it_cas2) THEN
759          frac=0.
760       ELSE
761          frac=(time_cas2-timeit)/(time_cas2-time_cas1)
762          frac=max(frac,0.0)
763       ENDIF
764
765       lat_prof_cas = lat_cas(it_cas2)                                   &
766     &          -frac*(lat_cas(it_cas2)-lat_cas(it_cas1))
767       sens_prof_cas = sens_cas(it_cas2)                                 &
768     &          -frac*(sens_cas(it_cas2)-sens_cas(it_cas1))
769       tkes_prof_cas = tkes_cas(it_cas2)                                   &
770     &          -frac*(tkes_cas(it_cas2)-tkes_cas(it_cas1))
771       ts_prof_cas = ts_cas(it_cas2)                                     &
772     &          -frac*(ts_cas(it_cas2)-ts_cas(it_cas1))
773       ps_prof_cas = ps_cas(it_cas2)                                     &
774     &          -frac*(ps_cas(it_cas2)-ps_cas(it_cas1))
775       ustar_prof_cas = ustar_cas(it_cas2)                               &
776     &          -frac*(ustar_cas(it_cas2)-ustar_cas(it_cas1))
777
778       do k=1,nlev_cas
779        plev_prof_cas(k) = plev_cas(k,it_cas2)                           &     
780     &          -frac*(plev_cas(k,it_cas2)-plev_cas(k,it_cas1))
781        t_prof_cas(k) = t_cas(k,it_cas2)                                 &       
782     &          -frac*(t_cas(k,it_cas2)-t_cas(k,it_cas1))
783        !print *,'k,frac,plev_cas1,plev_cas2=',k,frac,plev_cas(k,it_cas1),plev_cas(k,it_cas2)
784        theta_prof_cas(k) = theta_cas(k,it_cas2)                         &                     
785     &          -frac*(theta_cas(k,it_cas2)-theta_cas(k,it_cas1))
786        thv_prof_cas(k) = thv_cas(k,it_cas2)                             &         
787     &          -frac*(thv_cas(k,it_cas2)-thv_cas(k,it_cas1))
788        thl_prof_cas(k) = thl_cas(k,it_cas2)                             &             
789     &          -frac*(thl_cas(k,it_cas2)-thl_cas(k,it_cas1))
790        qv_prof_cas(k) = qv_cas(k,it_cas2)                               &
791     &          -frac*(qv_cas(k,it_cas2)-qv_cas(k,it_cas1))
792        ql_prof_cas(k) = ql_cas(k,it_cas2)                               &
793     &          -frac*(ql_cas(k,it_cas2)-ql_cas(k,it_cas1))
794        qi_prof_cas(k) = qi_cas(k,it_cas2)                               &
795     &          -frac*(qi_cas(k,it_cas2)-qi_cas(k,it_cas1))
796        u_prof_cas(k) = u_cas(k,it_cas2)                                 &
797     &          -frac*(u_cas(k,it_cas2)-u_cas(k,it_cas1))
798        v_prof_cas(k) = v_cas(k,it_cas2)                                 &
799     &          -frac*(v_cas(k,it_cas2)-v_cas(k,it_cas1))
800        ug_prof_cas(k) = ug_cas(k,it_cas2)                               &
801     &          -frac*(ug_cas(k,it_cas2)-ug_cas(k,it_cas1))
802        vg_prof_cas(k) = vg_cas(k,it_cas2)                               &
803     &          -frac*(vg_cas(k,it_cas2)-vg_cas(k,it_cas1))
804        temp_nudg_prof_cas(k) = temp_nudg_cas(k,it_cas2)                    &
805     &          -frac*(temp_nudg_cas(k,it_cas2)-temp_nudg_cas(k,it_cas1))
806        qv_nudg_prof_cas(k) = qv_nudg_cas(k,it_cas2)                        &
807     &          -frac*(qv_nudg_cas(k,it_cas2)-qv_nudg_cas(k,it_cas1))
808        u_nudg_prof_cas(k) = u_nudg_cas(k,it_cas2)                          &
809     &          -frac*(u_nudg_cas(k,it_cas2)-u_nudg_cas(k,it_cas1))
810        v_nudg_prof_cas(k) = v_nudg_cas(k,it_cas2)                          &
811     &          -frac*(v_nudg_cas(k,it_cas2)-v_nudg_cas(k,it_cas1))
812        vitw_prof_cas(k) = vitw_cas(k,it_cas2)                           &
813     &          -frac*(vitw_cas(k,it_cas2)-vitw_cas(k,it_cas1))
814        omega_prof_cas(k) = omega_cas(k,it_cas2)                         &
815     &          -frac*(omega_cas(k,it_cas2)-omega_cas(k,it_cas1))
816        tke_prof_cas(k) = tke_cas(k,it_cas2)                         &
817     &          -frac*(tke_cas(k,it_cas2)-tke_cas(k,it_cas1))
818        du_prof_cas(k) = du_cas(k,it_cas2)                               &
819     &          -frac*(du_cas(k,it_cas2)-du_cas(k,it_cas1))
820        hu_prof_cas(k) = hu_cas(k,it_cas2)                               &
821     &          -frac*(hu_cas(k,it_cas2)-hu_cas(k,it_cas1))
822        vu_prof_cas(k) = vu_cas(k,it_cas2)                               &
823     &          -frac*(vu_cas(k,it_cas2)-vu_cas(k,it_cas1))
824        dv_prof_cas(k) = dv_cas(k,it_cas2)                               &
825     &          -frac*(dv_cas(k,it_cas2)-dv_cas(k,it_cas1))
826        hv_prof_cas(k) = hv_cas(k,it_cas2)                               &
827     &          -frac*(hv_cas(k,it_cas2)-hv_cas(k,it_cas1))
828        vv_prof_cas(k) = vv_cas(k,it_cas2)                               &
829     &          -frac*(vv_cas(k,it_cas2)-vv_cas(k,it_cas1))
830        dt_prof_cas(k) = dt_cas(k,it_cas2)                               &
831     &          -frac*(dt_cas(k,it_cas2)-dt_cas(k,it_cas1))
832        ht_prof_cas(k) = ht_cas(k,it_cas2)                               &
833     &          -frac*(ht_cas(k,it_cas2)-ht_cas(k,it_cas1))
834        vt_prof_cas(k) = vt_cas(k,it_cas2)                               &
835     &          -frac*(vt_cas(k,it_cas2)-vt_cas(k,it_cas1))
836        dth_prof_cas(k) = dth_cas(k,it_cas2)                             &
837     &          -frac*(dth_cas(k,it_cas2)-dth_cas(k,it_cas1))
838        hth_prof_cas(k) = hth_cas(k,it_cas2)                             &
839     &          -frac*(hth_cas(k,it_cas2)-hth_cas(k,it_cas1))
840        vth_prof_cas(k) = vth_cas(k,it_cas2)                             &
841     &          -frac*(vth_cas(k,it_cas2)-vth_cas(k,it_cas1))
842        dtrad_prof_cas(k) = dtrad_cas(k,it_cas2)                         &
843     &          -frac*(dtrad_cas(k,it_cas2)-dtrad_cas(k,it_cas1))
844        dq_prof_cas(k) = dq_cas(k,it_cas2)                               &
845     &          -frac*(dq_cas(k,it_cas2)-dq_cas(k,it_cas1))
846        hq_prof_cas(k) = hq_cas(k,it_cas2)                               &
847     &          -frac*(hq_cas(k,it_cas2)-hq_cas(k,it_cas1))
848        vq_prof_cas(k) = vq_cas(k,it_cas2)                               &
849     &          -frac*(vq_cas(k,it_cas2)-vq_cas(k,it_cas1))
850       uw_prof_cas(k) = uw_cas(k,it_cas2)                                &
851     &          -frac*(uw_cas(k,it_cas2)-uw_cas(k,it_cas1))
852       vw_prof_cas(k) = vw_cas(k,it_cas2)                                &
853     &          -frac*(vw_cas(k,it_cas2)-vw_cas(k,it_cas1))
854       q1_prof_cas(k) = q1_cas(k,it_cas2)                                &
855     &          -frac*(q1_cas(k,it_cas2)-q1_cas(k,it_cas1))
856       q2_prof_cas(k) = q2_cas(k,it_cas2)                                &
857     &          -frac*(q2_cas(k,it_cas2)-q2_cas(k,it_cas1))
858        enddo
859
860        return
861        END SUBROUTINE interp_case_time_std
862
863!**********************************************************************************************
864!=====================================================================
865       SUBROUTINE interp2_case_vertical_std(play,plev,nlev_cas,plev_prof_cas                           &
866     &         ,t_prof_cas,th_prof_cas,thv_prof_cas,thl_prof_cas                                       &
867     &         ,qv_prof_cas,ql_prof_cas,qi_prof_cas,u_prof_cas,v_prof_cas                              &
868     &         ,ug_prof_cas,vg_prof_cas                                                                &
869     &         ,temp_nudg_prof_cas,qv_nudg_prof_cas,u_nudg_prof_cas,v_nudg_prof_cas                    &
870     &         ,vitw_prof_cas,omega_prof_cas,tke_prof_cas                                              &
871     &         ,du_prof_cas,hu_prof_cas,vu_prof_cas,dv_prof_cas,hv_prof_cas,vv_prof_cas                &
872     &         ,dt_prof_cas,ht_prof_cas,vt_prof_cas,dtrad_prof_cas,dq_prof_cas,hq_prof_cas,vq_prof_cas &
873     &         ,dth_prof_cas,hth_prof_cas,vth_prof_cas                                                 &
874!
875     &         ,t_mod_cas,theta_mod_cas,thv_mod_cas,thl_mod_cas                                        &
876     &         ,qv_mod_cas,ql_mod_cas,qi_mod_cas,u_mod_cas,v_mod_cas                                   &
877     &         ,ug_mod_cas,vg_mod_cas                                                                  &
878     &         ,temp_nudg_mod_cas,qv_nudg_mod_cas,u_nudg_mod_cas,v_nudg_mod_cas                        &
879     &         ,w_mod_cas,omega_mod_cas,tke_mod_cas                                                    &
880     &         ,du_mod_cas,hu_mod_cas,vu_mod_cas,dv_mod_cas,hv_mod_cas,vv_mod_cas                      &
881     &         ,dt_mod_cas,ht_mod_cas,vt_mod_cas,dtrad_mod_cas,dq_mod_cas,hq_mod_cas,vq_mod_cas        &
882     &         ,dth_mod_cas,hth_mod_cas,vth_mod_cas,mxcalc)
883 
884       implicit none
885 
886#include "YOMCST.h"
887#include "dimensions.h"
888
889!-------------------------------------------------------------------------
890! Vertical interpolation of generic case forcing data onto mod_casel levels
891!-------------------------------------------------------------------------
892 
893       integer nlevmax
894       parameter (nlevmax=41)
895       integer nlev_cas,mxcalc
896!       real play(llm), plev_prof(nlevmax)
897!       real t_prof(nlevmax),q_prof(nlevmax)
898!       real u_prof(nlevmax),v_prof(nlevmax), w_prof(nlevmax)
899!       real ht_prof(nlevmax),vt_prof(nlevmax)
900!       real hq_prof(nlevmax),vq_prof(nlevmax)
901 
902       real play(llm), plev(llm+1), plev_prof_cas(nlev_cas)
903       real t_prof_cas(nlev_cas),th_prof_cas(nlev_cas),thv_prof_cas(nlev_cas),thl_prof_cas(nlev_cas)
904       real qv_prof_cas(nlev_cas),ql_prof_cas(nlev_cas),qi_prof_cas(nlev_cas)
905       real u_prof_cas(nlev_cas),v_prof_cas(nlev_cas)
906       real ug_prof_cas(nlev_cas),vg_prof_cas(nlev_cas), vitw_prof_cas(nlev_cas),omega_prof_cas(nlev_cas),tke_prof_cas(nlev_cas)
907       real temp_nudg_prof_cas(nlev_cas),qv_nudg_prof_cas(nlev_cas)
908       real u_nudg_prof_cas(nlev_cas),v_nudg_prof_cas(nlev_cas)
909
910       real du_prof_cas(nlev_cas),hu_prof_cas(nlev_cas),vu_prof_cas(nlev_cas)
911       real dv_prof_cas(nlev_cas),hv_prof_cas(nlev_cas),vv_prof_cas(nlev_cas)
912       real dt_prof_cas(nlev_cas),ht_prof_cas(nlev_cas),vt_prof_cas(nlev_cas),dtrad_prof_cas(nlev_cas)
913       real dth_prof_cas(nlev_cas),hth_prof_cas(nlev_cas),vth_prof_cas(nlev_cas)
914       real dq_prof_cas(nlev_cas),hq_prof_cas(nlev_cas),vq_prof_cas(nlev_cas)
915 
916       real t_mod_cas(llm),theta_mod_cas(llm),thv_mod_cas(llm),thl_mod_cas(llm)
917       real qv_mod_cas(llm),ql_mod_cas(llm),qi_mod_cas(llm)
918       real u_mod_cas(llm),v_mod_cas(llm)
919       real ug_mod_cas(llm),vg_mod_cas(llm), w_mod_cas(llm),omega_mod_cas(llm),tke_mod_cas(llm+1)
920       real temp_nudg_mod_cas(llm),qv_nudg_mod_cas(llm)
921       real u_nudg_mod_cas(llm),v_nudg_mod_cas(llm)
922       real du_mod_cas(llm),hu_mod_cas(llm),vu_mod_cas(llm)
923       real dv_mod_cas(llm),hv_mod_cas(llm),vv_mod_cas(llm)
924       real dt_mod_cas(llm),ht_mod_cas(llm),vt_mod_cas(llm),dtrad_mod_cas(llm)
925       real dth_mod_cas(llm),hth_mod_cas(llm),vth_mod_cas(llm)
926       real dq_mod_cas(llm),hq_mod_cas(llm),vq_mod_cas(llm)
927 
928       integer l,k,k1,k2
929       real frac,frac1,frac2,fact
930 
931
932
933! for variables defined at the middle of layers
934
935       do l = 1, llm
936
937        if (play(l).ge.plev_prof_cas(nlev_cas)) then
938 
939        mxcalc=l
940!        print *,'debut interp2, mxcalc=',mxcalc
941         k1=0
942         k2=0
943
944         if (play(l).le.plev_prof_cas(1)) then
945
946         do k = 1, nlev_cas-1
947          if (play(l).le.plev_prof_cas(k).and. play(l).gt.plev_prof_cas(k+1)) then
948            k1=k
949            k2=k+1
950          endif
951         enddo
952
953         if (k1.eq.0 .or. k2.eq.0) then
954          write(*,*) 'PB! k1, k2 = ',k1,k2
955          write(*,*) 'l,play(l) = ',l,play(l)/100
956         do k = 1, nlev_cas-1
957          write(*,*) 'k,plev_prof_cas(k) = ',k,plev_prof_cas(k)/100
958         enddo
959         endif
960
961
962
963         frac = (plev_prof_cas(k2)-play(l))/(plev_prof_cas(k2)-plev_prof_cas(k1))
964         
965         t_mod_cas(l)= t_prof_cas(k2) - frac*(t_prof_cas(k2)-t_prof_cas(k1))
966         theta_mod_cas(l)= th_prof_cas(k2) - frac*(th_prof_cas(k2)-th_prof_cas(k1))
967         if(theta_mod_cas(l).NE.0) t_mod_cas(l)= theta_mod_cas(l)*(play(l)/100000.)**(RD/RCPD)
968         thv_mod_cas(l)= thv_prof_cas(k2) - frac*(thv_prof_cas(k2)-thv_prof_cas(k1))
969         thl_mod_cas(l)= thl_prof_cas(k2) - frac*(thl_prof_cas(k2)-thl_prof_cas(k1))
970         qv_mod_cas(l)= qv_prof_cas(k2) - frac*(qv_prof_cas(k2)-qv_prof_cas(k1))
971         ql_mod_cas(l)= ql_prof_cas(k2) - frac*(ql_prof_cas(k2)-ql_prof_cas(k1))
972         qi_mod_cas(l)= qi_prof_cas(k2) - frac*(qi_prof_cas(k2)-qi_prof_cas(k1))
973         u_mod_cas(l)= u_prof_cas(k2) - frac*(u_prof_cas(k2)-u_prof_cas(k1))
974         v_mod_cas(l)= v_prof_cas(k2) - frac*(v_prof_cas(k2)-v_prof_cas(k1))
975         ug_mod_cas(l)= ug_prof_cas(k2) - frac*(ug_prof_cas(k2)-ug_prof_cas(k1))
976         vg_mod_cas(l)= vg_prof_cas(k2) - frac*(vg_prof_cas(k2)-vg_prof_cas(k1))
977         temp_nudg_mod_cas(l)= temp_nudg_prof_cas(k2) - frac*(temp_nudg_prof_cas(k2)-temp_nudg_prof_cas(k1))
978         qv_nudg_mod_cas(l)= qv_nudg_prof_cas(k2) - frac*(qv_nudg_prof_cas(k2)-qv_nudg_prof_cas(k1))
979         u_nudg_mod_cas(l)= u_nudg_prof_cas(k2) - frac*(u_nudg_prof_cas(k2)-u_nudg_prof_cas(k1))
980         v_nudg_mod_cas(l)= v_nudg_prof_cas(k2) - frac*(v_nudg_prof_cas(k2)-v_nudg_prof_cas(k1))
981         w_mod_cas(l)= vitw_prof_cas(k2) - frac*(vitw_prof_cas(k2)-vitw_prof_cas(k1))
982         omega_mod_cas(l)= omega_prof_cas(k2) - frac*(omega_prof_cas(k2)-omega_prof_cas(k1))
983         du_mod_cas(l)= du_prof_cas(k2) - frac*(du_prof_cas(k2)-du_prof_cas(k1))
984         hu_mod_cas(l)= hu_prof_cas(k2) - frac*(hu_prof_cas(k2)-hu_prof_cas(k1))
985         vu_mod_cas(l)= vu_prof_cas(k2) - frac*(vu_prof_cas(k2)-vu_prof_cas(k1))
986         dv_mod_cas(l)= dv_prof_cas(k2) - frac*(dv_prof_cas(k2)-dv_prof_cas(k1))
987         hv_mod_cas(l)= hv_prof_cas(k2) - frac*(hv_prof_cas(k2)-hv_prof_cas(k1))
988         vv_mod_cas(l)= vv_prof_cas(k2) - frac*(vv_prof_cas(k2)-vv_prof_cas(k1))
989         dt_mod_cas(l)= dt_prof_cas(k2) - frac*(dt_prof_cas(k2)-dt_prof_cas(k1))
990         ht_mod_cas(l)= ht_prof_cas(k2) - frac*(ht_prof_cas(k2)-ht_prof_cas(k1))
991         vt_mod_cas(l)= vt_prof_cas(k2) - frac*(vt_prof_cas(k2)-vt_prof_cas(k1))
992         dth_mod_cas(l)= dth_prof_cas(k2) - frac*(dth_prof_cas(k2)-dth_prof_cas(k1))
993         hth_mod_cas(l)= hth_prof_cas(k2) - frac*(hth_prof_cas(k2)-hth_prof_cas(k1))
994         vth_mod_cas(l)= vth_prof_cas(k2) - frac*(vth_prof_cas(k2)-vth_prof_cas(k1))
995         dq_mod_cas(l)= dq_prof_cas(k2) - frac*(dq_prof_cas(k2)-dq_prof_cas(k1))
996         hq_mod_cas(l)= hq_prof_cas(k2) - frac*(hq_prof_cas(k2)-hq_prof_cas(k1))
997         vq_mod_cas(l)= vq_prof_cas(k2) - frac*(vq_prof_cas(k2)-vq_prof_cas(k1))
998         dtrad_mod_cas(l)= dtrad_prof_cas(k2) - frac*(dtrad_prof_cas(k2)-dtrad_prof_cas(k1))
999     
1000         else !play>plev_prof_cas(1)
1001
1002         k1=1
1003         k2=2
1004         print *,'interp2_vert, k1,k2=',plev_prof_cas(k1),plev_prof_cas(k2)
1005         frac1 = (play(l)-plev_prof_cas(k2))/(plev_prof_cas(k1)-plev_prof_cas(k2))
1006         frac2 = (play(l)-plev_prof_cas(k1))/(plev_prof_cas(k1)-plev_prof_cas(k2))
1007         t_mod_cas(l)= frac1*t_prof_cas(k1) - frac2*t_prof_cas(k2)
1008         theta_mod_cas(l)= frac1*th_prof_cas(k1) - frac2*th_prof_cas(k2)
1009         if(theta_mod_cas(l).NE.0) t_mod_cas(l)= theta_mod_cas(l)*(play(l)/100000.)**(RD/RCPD)
1010         thv_mod_cas(l)= frac1*thv_prof_cas(k1) - frac2*thv_prof_cas(k2)
1011         thl_mod_cas(l)= frac1*thl_prof_cas(k1) - frac2*thl_prof_cas(k2)
1012         qv_mod_cas(l)= frac1*qv_prof_cas(k1) - frac2*qv_prof_cas(k2)
1013         ql_mod_cas(l)= frac1*ql_prof_cas(k1) - frac2*ql_prof_cas(k2)
1014         qi_mod_cas(l)= frac1*qi_prof_cas(k1) - frac2*qi_prof_cas(k2)
1015         u_mod_cas(l)= frac1*u_prof_cas(k1) - frac2*u_prof_cas(k2)
1016         v_mod_cas(l)= frac1*v_prof_cas(k1) - frac2*v_prof_cas(k2)
1017         ug_mod_cas(l)= frac1*ug_prof_cas(k1) - frac2*ug_prof_cas(k2)
1018         vg_mod_cas(l)= frac1*vg_prof_cas(k1) - frac2*vg_prof_cas(k2)
1019         temp_nudg_mod_cas(l)= frac1*temp_nudg_prof_cas(k1) - frac2*temp_nudg_prof_cas(k2)
1020         qv_nudg_mod_cas(l)= frac1*qv_nudg_prof_cas(k1) - frac2*qv_nudg_prof_cas(k2)
1021         u_nudg_mod_cas(l)= frac1*u_nudg_prof_cas(k1) - frac2*u_nudg_prof_cas(k2)
1022         v_nudg_mod_cas(l)= frac1*v_nudg_prof_cas(k1) - frac2*v_nudg_prof_cas(k2)
1023         w_mod_cas(l)= frac1*vitw_prof_cas(k1) - frac2*vitw_prof_cas(k2)
1024         omega_mod_cas(l)= frac1*omega_prof_cas(k1) - frac2*omega_prof_cas(k2)
1025         du_mod_cas(l)= frac1*du_prof_cas(k1) - frac2*du_prof_cas(k2)
1026         hu_mod_cas(l)= frac1*hu_prof_cas(k1) - frac2*hu_prof_cas(k2)
1027         vu_mod_cas(l)= frac1*vu_prof_cas(k1) - frac2*vu_prof_cas(k2)
1028         dv_mod_cas(l)= frac1*dv_prof_cas(k1) - frac2*dv_prof_cas(k2)
1029         hv_mod_cas(l)= frac1*hv_prof_cas(k1) - frac2*hv_prof_cas(k2)
1030         vv_mod_cas(l)= frac1*vv_prof_cas(k1) - frac2*vv_prof_cas(k2)
1031         dt_mod_cas(l)= frac1*dt_prof_cas(k1) - frac2*dt_prof_cas(k2)
1032         ht_mod_cas(l)= frac1*ht_prof_cas(k1) - frac2*ht_prof_cas(k2)
1033         vt_mod_cas(l)= frac1*vt_prof_cas(k1) - frac2*vt_prof_cas(k2)
1034         dth_mod_cas(l)= frac1*dth_prof_cas(k1) - frac2*dth_prof_cas(k2)
1035         hth_mod_cas(l)= frac1*hth_prof_cas(k1) - frac2*hth_prof_cas(k2)
1036         vth_mod_cas(l)= frac1*vth_prof_cas(k1) - frac2*vth_prof_cas(k2)
1037         dq_mod_cas(l)= frac1*dq_prof_cas(k1) - frac2*dq_prof_cas(k2)
1038         hq_mod_cas(l)= frac1*hq_prof_cas(k1) - frac2*hq_prof_cas(k2)
1039         vq_mod_cas(l)= frac1*vq_prof_cas(k1) - frac2*vq_prof_cas(k2)
1040         dtrad_mod_cas(l)= frac1*dtrad_prof_cas(k1) - frac2*dtrad_prof_cas(k2)
1041
1042         endif ! play.le.plev_prof_cas(1)
1043
1044        else ! above max altitude of forcing file
1045 
1046!jyg
1047         fact=20.*(plev_prof_cas(nlev_cas)-play(l))/plev_prof_cas(nlev_cas) !jyg
1048         fact = max(fact,0.)                                           !jyg
1049         fact = exp(-fact)                                             !jyg
1050         t_mod_cas(l)= t_prof_cas(nlev_cas)                            !jyg
1051         theta_mod_cas(l)= th_prof_cas(nlev_cas)                       !jyg
1052         thv_mod_cas(l)= thv_prof_cas(nlev_cas)                        !jyg
1053         thl_mod_cas(l)= thl_prof_cas(nlev_cas)                        !jyg
1054         qv_mod_cas(l)= qv_prof_cas(nlev_cas)*fact                     !jyg
1055         ql_mod_cas(l)= ql_prof_cas(nlev_cas)*fact                     !jyg
1056         qi_mod_cas(l)= qi_prof_cas(nlev_cas)*fact                     !jyg
1057         u_mod_cas(l)= u_prof_cas(nlev_cas)*fact                       !jyg
1058         v_mod_cas(l)= v_prof_cas(nlev_cas)*fact                       !jyg
1059         ug_mod_cas(l)= ug_prof_cas(nlev_cas)                          !jyg
1060         vg_mod_cas(l)= vg_prof_cas(nlev_cas)                          !jyg
1061         temp_nudg_mod_cas(l)= temp_nudg_prof_cas(nlev_cas)            !jyg
1062         qv_nudg_mod_cas(l)= qv_nudg_prof_cas(nlev_cas)                !jyg
1063         u_nudg_mod_cas(l)= u_nudg_prof_cas(nlev_cas)                  !jyg
1064         v_nudg_mod_cas(l)= v_nudg_prof_cas(nlev_cas)                  !jyg
1065         thv_mod_cas(l)= thv_prof_cas(nlev_cas)                        !jyg
1066         w_mod_cas(l)= 0.0                                             !jyg
1067         omega_mod_cas(l)= 0.0                                         !jyg
1068         du_mod_cas(l)= du_prof_cas(nlev_cas)*fact
1069         hu_mod_cas(l)= hu_prof_cas(nlev_cas)*fact                     !jyg
1070         vu_mod_cas(l)= vu_prof_cas(nlev_cas)*fact                     !jyg
1071         dv_mod_cas(l)= dv_prof_cas(nlev_cas)*fact
1072         hv_mod_cas(l)= hv_prof_cas(nlev_cas)*fact                     !jyg
1073         vv_mod_cas(l)= vv_prof_cas(nlev_cas)*fact                     !jyg
1074         dt_mod_cas(l)= dt_prof_cas(nlev_cas)
1075         ht_mod_cas(l)= ht_prof_cas(nlev_cas)                          !jyg
1076         vt_mod_cas(l)= vt_prof_cas(nlev_cas)                          !jyg
1077         dth_mod_cas(l)= dth_prof_cas(nlev_cas)
1078         hth_mod_cas(l)= hth_prof_cas(nlev_cas)                        !jyg
1079         vth_mod_cas(l)= vth_prof_cas(nlev_cas)                        !jyg
1080         dq_mod_cas(l)= dq_prof_cas(nlev_cas)*fact
1081         hq_mod_cas(l)= hq_prof_cas(nlev_cas)*fact                     !jyg
1082         vq_mod_cas(l)= vq_prof_cas(nlev_cas)*fact                     !jyg
1083         dtrad_mod_cas(l)= dtrad_prof_cas(nlev_cas)*fact               !jyg
1084 
1085        endif ! play
1086 
1087       enddo ! l
1088
1089! for variables defined at layer interfaces (EV):
1090
1091
1092       do l = 1, llm+1
1093
1094        if (plev(l).ge.plev_prof_cas(nlev_cas)) then
1095
1096         mxcalc=l
1097         k1=0
1098         k2=0
1099
1100         if (plev(l).le.plev_prof_cas(1)) then
1101
1102         do k = 1, nlev_cas-1
1103          if (plev(l).le.plev_prof_cas(k).and. plev(l).gt.plev_prof_cas(k+1)) then
1104            k1=k
1105            k2=k+1
1106          endif
1107         enddo
1108
1109         if (k1.eq.0 .or. k2.eq.0) then
1110          write(*,*) 'PB! k1, k2 = ',k1,k2
1111          write(*,*) 'l,plev(l) = ',l,plev(l)/100
1112         do k = 1, nlev_cas-1
1113          write(*,*) 'k,plev_prof_cas(k) = ',k,plev_prof_cas(k)/100
1114         enddo
1115         endif
1116
1117         frac = (plev_prof_cas(k2)-plev(l))/(plev_prof_cas(k2)-plev_prof_cas(k1))
1118         tke_mod_cas(l)= tke_prof_cas(k2) - frac*(tke_prof_cas(k2)-tke_prof_cas(k1))
1119         else !play>plev_prof_cas(1)
1120         k1=1
1121         k2=2
1122         tke_mod_cas(l)= frac1*tke_prof_cas(k1) - frac2*tke_prof_cas(k2)
1123
1124         endif ! plev.le.plev_prof_cas(1)
1125
1126        else ! above max altitude of forcing file
1127
1128         tke_mod_cas(l)=0.0
1129
1130        endif ! plev
1131
1132       enddo ! l
1133
1134
1135
1136          return
1137          end SUBROUTINE interp2_case_vertical_std
1138!*****************************************************************************
1139
1140
1141
1142
1143
1144END MODULE mod_1D_cases_read_std
Note: See TracBrowser for help on using the repository browser.