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

Last change on this file since 3669 was 3669, checked in by fhourdin, 5 years ago

Traitement du guidage pour le format standard
Marie-Pierre et Frédéric

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