source: LMDZ6/trunk/libf/phylmd/dyn1d/mod_1D_cases_read.F90 @ 4257

Last change on this file since 4257 was 2373, checked in by jyg, 9 years ago

From MPL:
Correction of some missing initializations in
dyn1D and of two 1D case bugs.

  • Property svn:keywords set to Id
File size: 38.0 KB
Line 
1!
2! $Id: mod_1D_cases_read.F90 2373 2015-10-13 17:28:01Z lguez $
3!
4MODULE mod_1D_cases_read
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!       integer year_ini_cas, day_ini_cas, mth_ini_cas
14!       real heure_ini_cas
15!       real day_ju_ini_cas   ! Julian day of case first day
16!       parameter (year_ini_cas=2011)
17!       parameter (year_ini_cas=1969)
18!       parameter (mth_ini_cas=10)
19!       parameter (mth_ini_cas=6)
20!       parameter (day_ini_cas=1)  ! 10 = 10Juil2006
21!       parameter (day_ini_cas=24)  ! 24 = 24 juin 1969
22!       parameter (heure_ini_cas=0.) !0h en secondes
23!       real pdt_cas
24!       parameter (pdt_cas=3.*3600)
25
26!CR ATTENTION TEST AMMA
27!        parameter (year_ini_cas=2006)
28!        parameter (mth_ini_cas=7)
29!        parameter (day_ini_cas=10)  ! 10 = 10Juil2006
30!        parameter (heure_ini_cas=0.) !0h en secondes
31!        parameter (pdt_cas=1800.)
32
33!profils environnementaux
34        real, allocatable::  plev_cas(:,:)
35
36        real, allocatable::  z_cas(:,:)
37        real, allocatable::  t_cas(:,:),q_cas(:,:),rh_cas(:,:)
38        real, allocatable::  th_cas(:,:),rv_cas(:,:)
39        real, allocatable::  u_cas(:,:)
40        real, allocatable::  v_cas(:,:)
41
42!forcing
43        real, allocatable::  ht_cas(:,:),vt_cas(:,:),dt_cas(:,:),dtrad_cas(:,:)
44        real, allocatable::  hth_cas(:,:),vth_cas(:,:),dth_cas(:,:)
45        real, allocatable::  hq_cas(:,:),vq_cas(:,:),dq_cas(:,:)
46        real, allocatable::  hr_cas(:,:),vr_cas(:,:),dr_cas(:,:)
47        real, allocatable::  hu_cas(:,:),vu_cas(:,:),du_cas(:,:)
48        real, allocatable::  hv_cas(:,:),vv_cas(:,:),dv_cas(:,:)
49        real, allocatable::  vitw_cas(:,:)
50        real, allocatable::  ug_cas(:,:),vg_cas(:,:)
51        real, allocatable::  lat_cas(:),sens_cas(:),ts_cas(:),ustar_cas(:)
52        real, allocatable::  uw_cas(:,:),vw_cas(:,:),q1_cas(:,:),q2_cas(:,:)
53
54!champs interpoles
55        real, allocatable::  plev_prof_cas(:)
56        real, allocatable::  t_prof_cas(:)
57        real, allocatable::  q_prof_cas(:)
58        real, allocatable::  u_prof_cas(:)
59        real, allocatable::  v_prof_cas(:)       
60
61        real, allocatable::  vitw_prof_cas(:)
62        real, allocatable::  ug_prof_cas(:)
63        real, allocatable::  vg_prof_cas(:)
64        real, allocatable::  ht_prof_cas(:)
65        real, allocatable::  hq_prof_cas(:)
66        real, allocatable::  vt_prof_cas(:)
67        real, allocatable::  vq_prof_cas(:)
68        real, allocatable::  dt_prof_cas(:)
69        real, allocatable::  dtrad_prof_cas(:)
70        real, allocatable::  dq_prof_cas(:)
71        real, allocatable::  hu_prof_cas(:)
72        real, allocatable::  hv_prof_cas(:)
73        real, allocatable::  vu_prof_cas(:)
74        real, allocatable::  vv_prof_cas(:)
75        real, allocatable::  du_prof_cas(:)
76        real, allocatable::  dv_prof_cas(:)
77        real, allocatable::  uw_prof_cas(:)
78        real, allocatable::  vw_prof_cas(:)
79        real, allocatable::  q1_prof_cas(:)
80        real, allocatable::  q2_prof_cas(:)
81
82
83        real lat_prof_cas,sens_prof_cas,ts_prof_cas,ustar_prof_cas
84     
85
86
87CONTAINS
88
89SUBROUTINE read_1D_cas
90      implicit none
91
92#include "netcdf.inc"
93
94      INTEGER nid,rid,ierr
95      INTEGER ii,jj
96
97      fich_cas='setup/cas.nc'
98      print*,'fich_cas ',fich_cas
99      ierr = NF_OPEN(fich_cas,NF_NOWRITE,nid)
100      print*,'fich_cas,NF_NOWRITE,nid ',fich_cas,NF_NOWRITE,nid
101      if (ierr.NE.NF_NOERR) then
102         write(*,*) 'ERROR: GROS Pb opening forcings nc file '
103         write(*,*) NF_STRERROR(ierr)
104         stop ""
105      endif
106!.......................................................................
107      ierr=NF_INQ_DIMID(nid,'lat',rid)
108      IF (ierr.NE.NF_NOERR) THEN
109         print*, 'Oh probleme lecture dimension lat'
110      ENDIF
111      ierr=NF_INQ_DIMLEN(nid,rid,ii)
112      print*,'OK1 nid,rid,lat',nid,rid,ii
113!.......................................................................
114      ierr=NF_INQ_DIMID(nid,'lon',rid)
115      IF (ierr.NE.NF_NOERR) THEN
116         print*, 'Oh probleme lecture dimension lon'
117      ENDIF
118      ierr=NF_INQ_DIMLEN(nid,rid,jj)
119      print*,'OK2 nid,rid,lat',nid,rid,jj
120!.......................................................................
121      ierr=NF_INQ_DIMID(nid,'lev',rid)
122      IF (ierr.NE.NF_NOERR) THEN
123         print*, 'Oh probleme lecture dimension zz'
124      ENDIF
125      ierr=NF_INQ_DIMLEN(nid,rid,nlev_cas)
126      print*,'OK3 nid,rid,nlev_cas',nid,rid,nlev_cas
127!.......................................................................
128      ierr=NF_INQ_DIMID(nid,'time',rid)
129      print*,'nid,rid',nid,rid
130      nt_cas=0
131      IF (ierr.NE.NF_NOERR) THEN
132        stop 'probleme lecture dimension sens'
133      ENDIF
134      ierr=NF_INQ_DIMLEN(nid,rid,nt_cas)
135      print*,'OK4 nid,rid,nt_cas',nid,rid,nt_cas
136
137!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
138!profils moyens:
139        allocate(plev_cas(nlev_cas,nt_cas))       
140        allocate(z_cas(nlev_cas,nt_cas))
141        allocate(t_cas(nlev_cas,nt_cas),q_cas(nlev_cas,nt_cas),rh_cas(nlev_cas,nt_cas))
142        allocate(th_cas(nlev_cas,nt_cas),rv_cas(nlev_cas,nt_cas))
143        allocate(u_cas(nlev_cas,nt_cas))
144        allocate(v_cas(nlev_cas,nt_cas))
145
146!forcing
147        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))
148        allocate(hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas),dq_cas(nlev_cas,nt_cas))
149        allocate(hth_cas(nlev_cas,nt_cas),vth_cas(nlev_cas,nt_cas),dth_cas(nlev_cas,nt_cas))
150        allocate(hr_cas(nlev_cas,nt_cas),vr_cas(nlev_cas,nt_cas),dr_cas(nlev_cas,nt_cas))
151        allocate(hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas),du_cas(nlev_cas,nt_cas))
152        allocate(hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas),dv_cas(nlev_cas,nt_cas))
153        allocate(vitw_cas(nlev_cas,nt_cas))
154        allocate(ug_cas(nlev_cas,nt_cas))
155        allocate(vg_cas(nlev_cas,nt_cas))
156        allocate(lat_cas(nt_cas),sens_cas(nt_cas),ts_cas(nt_cas),ustar_cas(nt_cas))
157        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))
158
159
160!champs interpoles
161        allocate(plev_prof_cas(nlev_cas))
162        allocate(t_prof_cas(nlev_cas))
163        allocate(q_prof_cas(nlev_cas))
164        allocate(u_prof_cas(nlev_cas))
165        allocate(v_prof_cas(nlev_cas))
166
167        allocate(vitw_prof_cas(nlev_cas))
168        allocate(ug_prof_cas(nlev_cas))
169        allocate(vg_prof_cas(nlev_cas))
170        allocate(ht_prof_cas(nlev_cas))
171        allocate(hq_prof_cas(nlev_cas))
172        allocate(hu_prof_cas(nlev_cas))
173        allocate(hv_prof_cas(nlev_cas))
174        allocate(vt_prof_cas(nlev_cas))
175        allocate(vq_prof_cas(nlev_cas))
176        allocate(vu_prof_cas(nlev_cas))
177        allocate(vv_prof_cas(nlev_cas))
178        allocate(dt_prof_cas(nlev_cas))
179        allocate(dtrad_prof_cas(nlev_cas))
180        allocate(dq_prof_cas(nlev_cas))
181        allocate(du_prof_cas(nlev_cas))
182        allocate(dv_prof_cas(nlev_cas))
183        allocate(uw_prof_cas(nlev_cas))
184        allocate(vw_prof_cas(nlev_cas))
185        allocate(q1_prof_cas(nlev_cas))
186        allocate(q2_prof_cas(nlev_cas))
187
188        print*,'Allocations OK'
189        call read_cas(nid,nlev_cas,nt_cas                                       &
190     &     ,z_cas,plev_cas,t_cas,q_cas,rh_cas,th_cas,rv_cas,u_cas,v_cas         &
191     &     ,ug_cas,vg_cas,vitw_cas,du_cas,hu_cas,vu_cas,dv_cas,hv_cas,vv_cas    &
192     &     ,dt_cas,dtrad_cas,ht_cas,vt_cas,dq_cas,hq_cas,vq_cas                 &
193     &     ,dth_cas,hth_cas,vth_cas,dr_cas,hr_cas,vr_cas,sens_cas,lat_cas,ts_cas&
194     &     ,ustar_cas,uw_cas,vw_cas,q1_cas,q2_cas)
195        print*,'Read cas OK'
196
197
198END SUBROUTINE read_1D_cas
199
200
201
202!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
203SUBROUTINE deallocate_1D_cases
204!profils environnementaux:
205        deallocate(plev_cas)
206       
207        deallocate(z_cas)
208        deallocate(t_cas,q_cas,rh_cas)
209        deallocate(th_cas,rv_cas)
210        deallocate(u_cas)
211        deallocate(v_cas)
212       
213!forcing
214        deallocate(ht_cas,vt_cas,dt_cas,dtrad_cas)
215        deallocate(hq_cas,vq_cas,dq_cas)
216        deallocate(hth_cas,vth_cas,dth_cas)
217        deallocate(hr_cas,vr_cas,dr_cas)
218        deallocate(hu_cas,vu_cas,du_cas)
219        deallocate(hv_cas,vv_cas,dv_cas)
220        deallocate(vitw_cas)
221        deallocate(ug_cas)
222        deallocate(vg_cas)
223        deallocate(lat_cas,sens_cas,ts_cas,ustar_cas,uw_cas,vw_cas,q1_cas,q2_cas)
224
225!champs interpoles
226        deallocate(plev_prof_cas)
227        deallocate(t_prof_cas)
228        deallocate(q_prof_cas)
229        deallocate(u_prof_cas)
230        deallocate(v_prof_cas)
231
232        deallocate(vitw_prof_cas)
233        deallocate(ug_prof_cas)
234        deallocate(vg_prof_cas)
235        deallocate(ht_prof_cas)
236        deallocate(hq_prof_cas)
237        deallocate(hu_prof_cas)
238        deallocate(hv_prof_cas)
239        deallocate(vt_prof_cas)
240        deallocate(vq_prof_cas)
241        deallocate(vu_prof_cas)
242        deallocate(vv_prof_cas)
243        deallocate(dt_prof_cas)
244        deallocate(dtrad_prof_cas)
245        deallocate(dq_prof_cas)
246        deallocate(du_prof_cas)
247        deallocate(dv_prof_cas)
248        deallocate(t_prof_cas)
249        deallocate(q_prof_cas)
250        deallocate(u_prof_cas)
251        deallocate(v_prof_cas)
252        deallocate(uw_prof_cas)
253        deallocate(vw_prof_cas)
254        deallocate(q1_prof_cas)
255        deallocate(q2_prof_cas)
256
257END SUBROUTINE deallocate_1D_cases
258
259
260END MODULE mod_1D_cases_read
261!=====================================================================
262      subroutine read_cas(nid,nlevel,ntime                          &
263     &     ,zz,pp,temp,qv,rh,theta,rv,u,v,ug,vg,w,                   &
264     &     du,hu,vu,dv,hv,vv,dt,dtrad,ht,vt,dq,hq,vq,                     &
265     &     dth,hth,vth,dr,hr,vr,sens,flat,ts,ustar,uw,vw,q1,q2)
266
267!program reading forcing of the case study
268      implicit none
269#include "netcdf.inc"
270
271      integer ntime,nlevel
272
273      real zz(nlevel,ntime)
274      real pp(nlevel,ntime)
275      real temp(nlevel,ntime),qv(nlevel,ntime),rh(nlevel,ntime)
276      real theta(nlevel,ntime),rv(nlevel,ntime)
277      real u(nlevel,ntime)
278      real v(nlevel,ntime)
279      real ug(nlevel,ntime)
280      real vg(nlevel,ntime)
281      real w(nlevel,ntime)
282      real du(nlevel,ntime),hu(nlevel,ntime),vu(nlevel,ntime)
283      real dv(nlevel,ntime),hv(nlevel,ntime),vv(nlevel,ntime)
284      real dt(nlevel,ntime),ht(nlevel,ntime),vt(nlevel,ntime)
285      real dtrad(nlevel,ntime)
286      real dq(nlevel,ntime),hq(nlevel,ntime),vq(nlevel,ntime)
287      real dth(nlevel,ntime),hth(nlevel,ntime),vth(nlevel,ntime)
288      real dr(nlevel,ntime),hr(nlevel,ntime),vr(nlevel,ntime)
289      real flat(ntime),sens(ntime),ts(ntime),ustar(ntime)
290      real uw(nlevel,ntime),vw(nlevel,ntime),q1(nlevel,ntime),q2(nlevel,ntime)
291
292
293      integer nid, ierr,rid
294      integer nbvar3d
295      parameter(nbvar3d=39)
296      integer var3didin(nbvar3d)
297
298       ierr=NF_INQ_VARID(nid,"zz",var3didin(1))
299         if(ierr/=NF_NOERR) then
300           write(*,*) NF_STRERROR(ierr)
301           stop 'lev'
302         endif
303     
304      ierr=NF_INQ_VARID(nid,"pp",var3didin(2))
305         if(ierr/=NF_NOERR) then
306           write(*,*) NF_STRERROR(ierr)
307           stop 'plev'
308         endif
309
310
311      ierr=NF_INQ_VARID(nid,"temp",var3didin(3))
312         if(ierr/=NF_NOERR) then
313           write(*,*) NF_STRERROR(ierr)
314           stop 'temp'
315         endif
316
317      ierr=NF_INQ_VARID(nid,"qv",var3didin(4))
318         if(ierr/=NF_NOERR) then
319           write(*,*) NF_STRERROR(ierr)
320           stop 'qv'
321         endif
322
323      ierr=NF_INQ_VARID(nid,"rh",var3didin(5))
324         if(ierr/=NF_NOERR) then
325           write(*,*) NF_STRERROR(ierr)
326           stop 'rh'
327         endif
328
329      ierr=NF_INQ_VARID(nid,"theta",var3didin(6))
330         if(ierr/=NF_NOERR) then
331           write(*,*) NF_STRERROR(ierr)
332           stop 'theta'
333         endif
334
335      ierr=NF_INQ_VARID(nid,"rv",var3didin(7))
336         if(ierr/=NF_NOERR) then
337           write(*,*) NF_STRERROR(ierr)
338           stop 'rv'
339         endif
340
341
342      ierr=NF_INQ_VARID(nid,"u",var3didin(8))
343         if(ierr/=NF_NOERR) then
344           write(*,*) NF_STRERROR(ierr)
345           stop 'u'
346         endif
347
348      ierr=NF_INQ_VARID(nid,"v",var3didin(9))
349         if(ierr/=NF_NOERR) then
350           write(*,*) NF_STRERROR(ierr)
351           stop 'v'
352         endif
353
354       ierr=NF_INQ_VARID(nid,"ug",var3didin(10))
355         if(ierr/=NF_NOERR) then
356           write(*,*) NF_STRERROR(ierr)
357           stop 'ug'
358         endif
359
360      ierr=NF_INQ_VARID(nid,"vg",var3didin(11))
361         if(ierr/=NF_NOERR) then
362           write(*,*) NF_STRERROR(ierr)
363           stop 'vg'
364         endif
365
366      ierr=NF_INQ_VARID(nid,"w",var3didin(12))
367         if(ierr/=NF_NOERR) then
368           write(*,*) NF_STRERROR(ierr)
369           stop 'w'
370         endif
371
372      ierr=NF_INQ_VARID(nid,"advu",var3didin(13))
373         if(ierr/=NF_NOERR) then
374           write(*,*) NF_STRERROR(ierr)
375           stop 'advu'
376         endif
377
378      ierr=NF_INQ_VARID(nid,"hu",var3didin(14))
379         if(ierr/=NF_NOERR) then
380           write(*,*) NF_STRERROR(ierr)
381           stop 'hu'
382         endif
383
384       ierr=NF_INQ_VARID(nid,"vu",var3didin(15))
385         if(ierr/=NF_NOERR) then
386           write(*,*) NF_STRERROR(ierr)
387           stop 'vu'
388         endif
389
390       ierr=NF_INQ_VARID(nid,"advv",var3didin(16))
391         if(ierr/=NF_NOERR) then
392           write(*,*) NF_STRERROR(ierr)
393           stop 'advv'
394         endif
395
396      ierr=NF_INQ_VARID(nid,"hv",var3didin(17))
397         if(ierr/=NF_NOERR) then
398           write(*,*) NF_STRERROR(ierr)
399           stop 'hv'
400         endif
401
402       ierr=NF_INQ_VARID(nid,"vv",var3didin(18))
403         if(ierr/=NF_NOERR) then
404           write(*,*) NF_STRERROR(ierr)
405           stop 'vv'
406         endif
407
408      ierr=NF_INQ_VARID(nid,"advT",var3didin(19))
409         if(ierr/=NF_NOERR) then
410           write(*,*) NF_STRERROR(ierr)
411           stop 'advT'
412         endif
413
414      ierr=NF_INQ_VARID(nid,"hT",var3didin(20))
415         if(ierr/=NF_NOERR) then
416           write(*,*) NF_STRERROR(ierr)
417           stop 'hT'
418         endif
419
420      ierr=NF_INQ_VARID(nid,"vT",var3didin(21))
421         if(ierr/=NF_NOERR) then
422           write(*,*) NF_STRERROR(ierr)
423           stop 'vT'
424         endif
425
426      ierr=NF_INQ_VARID(nid,"advq",var3didin(22))
427         if(ierr/=NF_NOERR) then
428           write(*,*) NF_STRERROR(ierr)
429           stop 'advq'
430         endif
431     
432      ierr=NF_INQ_VARID(nid,"hq",var3didin(23))
433         if(ierr/=NF_NOERR) then
434           write(*,*) NF_STRERROR(ierr)
435           stop 'hq'
436         endif
437
438      ierr=NF_INQ_VARID(nid,"vq",var3didin(24))
439         if(ierr/=NF_NOERR) then
440           write(*,*) NF_STRERROR(ierr)
441           stop 'vq'
442         endif
443
444      ierr=NF_INQ_VARID(nid,"advth",var3didin(25))
445         if(ierr/=NF_NOERR) then
446           write(*,*) NF_STRERROR(ierr)
447           stop 'advth'
448         endif
449
450      ierr=NF_INQ_VARID(nid,"hth",var3didin(26))
451         if(ierr/=NF_NOERR) then
452           write(*,*) NF_STRERROR(ierr)
453           stop 'hth'
454         endif
455
456      ierr=NF_INQ_VARID(nid,"vth",var3didin(27))
457         if(ierr/=NF_NOERR) then
458           write(*,*) NF_STRERROR(ierr)
459           stop 'vth'
460         endif
461
462      ierr=NF_INQ_VARID(nid,"advr",var3didin(28))
463         if(ierr/=NF_NOERR) then
464           write(*,*) NF_STRERROR(ierr)
465           stop 'advr'
466         endif
467     
468      ierr=NF_INQ_VARID(nid,"hr",var3didin(29))
469         if(ierr/=NF_NOERR) then
470           write(*,*) NF_STRERROR(ierr)
471           stop 'hr'
472         endif
473
474      ierr=NF_INQ_VARID(nid,"vr",var3didin(30))
475         if(ierr/=NF_NOERR) then
476           write(*,*) NF_STRERROR(ierr)
477           stop 'vr'
478         endif
479
480      ierr=NF_INQ_VARID(nid,"radT",var3didin(31))
481         if(ierr/=NF_NOERR) then
482           write(*,*) NF_STRERROR(ierr)
483           stop 'radT'
484         endif
485
486      ierr=NF_INQ_VARID(nid,"sens",var3didin(32))
487         if(ierr/=NF_NOERR) then
488           write(*,*) NF_STRERROR(ierr)
489           stop 'sens'
490         endif
491
492      ierr=NF_INQ_VARID(nid,"flat",var3didin(33))
493         if(ierr/=NF_NOERR) then
494           write(*,*) NF_STRERROR(ierr)
495           stop 'flat'
496         endif
497
498      ierr=NF_INQ_VARID(nid,"ts",var3didin(34))
499         if(ierr/=NF_NOERR) then
500           write(*,*) NF_STRERROR(ierr)
501           stop 'ts'
502         endif
503
504      ierr=NF_INQ_VARID(nid,"ustar",var3didin(35))
505         if(ierr/=NF_NOERR) then
506           write(*,*) NF_STRERROR(ierr)
507           stop 'ustar'
508         endif
509
510      ierr=NF_INQ_VARID(nid,"uw",var3didin(36))
511         if(ierr/=NF_NOERR) then
512           write(*,*) NF_STRERROR(ierr)
513           stop 'uw'
514         endif
515
516      ierr=NF_INQ_VARID(nid,"vw",var3didin(37))
517         if(ierr/=NF_NOERR) then
518           write(*,*) NF_STRERROR(ierr)
519           stop 'vw'
520         endif
521
522      ierr=NF_INQ_VARID(nid,"q1",var3didin(38))
523         if(ierr/=NF_NOERR) then
524           write(*,*) NF_STRERROR(ierr)
525           stop 'q1'
526         endif
527
528      ierr=NF_INQ_VARID(nid,"q2",var3didin(39))
529         if(ierr/=NF_NOERR) then
530           write(*,*) NF_STRERROR(ierr)
531           stop 'q2'
532         endif
533 
534#ifdef NC_DOUBLE
535         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(1),zz)
536#else
537         ierr = NF_GET_VAR_REAL(nid,var3didin(1),zz)
538#endif
539         if(ierr/=NF_NOERR) then
540            write(*,*) NF_STRERROR(ierr)
541            stop "getvarup"
542         endif
543!          write(*,*)'lecture z ok',zz
544
545#ifdef NC_DOUBLE
546         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(2),pp)
547#else
548         ierr = NF_GET_VAR_REAL(nid,var3didin(2),pp)
549#endif
550         if(ierr/=NF_NOERR) then
551            write(*,*) NF_STRERROR(ierr)
552            stop "getvarup"
553         endif
554!          write(*,*)'lecture pp ok',pp
555
556
557#ifdef NC_DOUBLE
558         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(3),temp)
559#else
560         ierr = NF_GET_VAR_REAL(nid,var3didin(3),temp)
561#endif
562         if(ierr/=NF_NOERR) then
563            write(*,*) NF_STRERROR(ierr)
564            stop "getvarup"
565         endif
566!          write(*,*)'lecture T ok',temp
567
568#ifdef NC_DOUBLE
569         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(4),qv)
570#else
571         ierr = NF_GET_VAR_REAL(nid,var3didin(4),qv)
572#endif
573         if(ierr/=NF_NOERR) then
574            write(*,*) NF_STRERROR(ierr)
575            stop "getvarup"
576         endif
577!          write(*,*)'lecture qv ok',qv
578 
579#ifdef NC_DOUBLE
580         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(5),rh)
581#else
582         ierr = NF_GET_VAR_REAL(nid,var3didin(5),rh)
583#endif
584         if(ierr/=NF_NOERR) then
585            write(*,*) NF_STRERROR(ierr)
586            stop "getvarup"
587         endif
588!          write(*,*)'lecture rh ok',rh
589
590#ifdef NC_DOUBLE
591         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(6),theta)
592#else
593         ierr = NF_GET_VAR_REAL(nid,var3didin(6),theta)
594#endif
595         if(ierr/=NF_NOERR) then
596            write(*,*) NF_STRERROR(ierr)
597            stop "getvarup"
598         endif
599!          write(*,*)'lecture theta ok',theta
600
601#ifdef NC_DOUBLE
602         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(7),rv)
603#else
604         ierr = NF_GET_VAR_REAL(nid,var3didin(7),rv)
605#endif
606         if(ierr/=NF_NOERR) then
607            write(*,*) NF_STRERROR(ierr)
608            stop "getvarup"
609         endif
610!          write(*,*)'lecture rv ok',rv
611
612#ifdef NC_DOUBLE
613         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(8),u)
614#else
615         ierr = NF_GET_VAR_REAL(nid,var3didin(8),u)
616#endif
617         if(ierr/=NF_NOERR) then
618            write(*,*) NF_STRERROR(ierr)
619            stop "getvarup"
620         endif
621!          write(*,*)'lecture u ok',u
622
623#ifdef NC_DOUBLE
624         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(9),v)
625#else
626         ierr = NF_GET_VAR_REAL(nid,var3didin(9),v)
627#endif
628         if(ierr/=NF_NOERR) then
629            write(*,*) NF_STRERROR(ierr)
630            stop "getvarup"
631         endif
632!          write(*,*)'lecture v ok',v
633
634#ifdef NC_DOUBLE
635         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(10),ug)
636#else
637         ierr = NF_GET_VAR_REAL(nid,var3didin(10),ug)
638#endif
639         if(ierr/=NF_NOERR) then
640            write(*,*) NF_STRERROR(ierr)
641            stop "getvarup"
642         endif
643!          write(*,*)'lecture ug ok',ug
644
645#ifdef NC_DOUBLE
646         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(11),vg)
647#else
648         ierr = NF_GET_VAR_REAL(nid,var3didin(11),vg)
649#endif
650         if(ierr/=NF_NOERR) then
651            write(*,*) NF_STRERROR(ierr)
652            stop "getvarup"
653         endif
654!          write(*,*)'lecture vg ok',vg
655
656#ifdef NC_DOUBLE
657         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(12),w)
658#else
659         ierr = NF_GET_VAR_REAL(nid,var3didin(12),w)
660#endif
661         if(ierr/=NF_NOERR) then
662            write(*,*) NF_STRERROR(ierr)
663            stop "getvarup"
664         endif
665!          write(*,*)'lecture w ok',w
666
667#ifdef NC_DOUBLE
668         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(13),du)
669#else
670         ierr = NF_GET_VAR_REAL(nid,var3didin(13),du)
671#endif
672         if(ierr/=NF_NOERR) then
673            write(*,*) NF_STRERROR(ierr)
674            stop "getvarup"
675         endif
676!          write(*,*)'lecture du ok',du
677
678#ifdef NC_DOUBLE
679         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(14),hu)
680#else
681         ierr = NF_GET_VAR_REAL(nid,var3didin(14),hu)
682#endif
683         if(ierr/=NF_NOERR) then
684            write(*,*) NF_STRERROR(ierr)
685            stop "getvarup"
686         endif
687!          write(*,*)'lecture hu ok',hu
688
689#ifdef NC_DOUBLE
690         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(15),vu)
691#else
692         ierr = NF_GET_VAR_REAL(nid,var3didin(15),vu)
693#endif
694         if(ierr/=NF_NOERR) then
695            write(*,*) NF_STRERROR(ierr)
696            stop "getvarup"
697         endif
698!          write(*,*)'lecture vu ok',vu
699
700#ifdef NC_DOUBLE
701         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(16),dv)
702#else
703         ierr = NF_GET_VAR_REAL(nid,var3didin(16),dv)
704#endif
705         if(ierr/=NF_NOERR) then
706            write(*,*) NF_STRERROR(ierr)
707            stop "getvarup"
708         endif
709!          write(*,*)'lecture dv ok',dv
710
711#ifdef NC_DOUBLE
712         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(17),hv)
713#else
714         ierr = NF_GET_VAR_REAL(nid,var3didin(17),hv)
715#endif
716         if(ierr/=NF_NOERR) then
717            write(*,*) NF_STRERROR(ierr)
718            stop "getvarup"
719         endif
720!          write(*,*)'lecture hv ok',hv
721
722#ifdef NC_DOUBLE
723         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(18),vv)
724#else
725         ierr = NF_GET_VAR_REAL(nid,var3didin(18),vv)
726#endif
727         if(ierr/=NF_NOERR) then
728            write(*,*) NF_STRERROR(ierr)
729            stop "getvarup"
730         endif
731!          write(*,*)'lecture vv ok',vv
732
733#ifdef NC_DOUBLE
734         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(19),dt)
735#else
736         ierr = NF_GET_VAR_REAL(nid,var3didin(19),dt)
737#endif
738         if(ierr/=NF_NOERR) then
739            write(*,*) NF_STRERROR(ierr)
740            stop "getvarup"
741         endif
742!          write(*,*)'lecture dt ok',dt
743
744#ifdef NC_DOUBLE
745         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(20),ht)
746#else
747         ierr = NF_GET_VAR_REAL(nid,var3didin(20),ht)
748#endif
749         if(ierr/=NF_NOERR) then
750            write(*,*) NF_STRERROR(ierr)
751            stop "getvarup"
752         endif
753!          write(*,*)'lecture ht ok',ht
754
755#ifdef NC_DOUBLE
756         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(21),vt)
757#else
758         ierr = NF_GET_VAR_REAL(nid,var3didin(21),vt)
759#endif
760         if(ierr/=NF_NOERR) then
761            write(*,*) NF_STRERROR(ierr)
762            stop "getvarup"
763         endif
764!          write(*,*)'lecture vt ok',vt
765
766#ifdef NC_DOUBLE
767         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(22),dq)
768#else
769         ierr = NF_GET_VAR_REAL(nid,var3didin(22),dq)
770#endif
771         if(ierr/=NF_NOERR) then
772            write(*,*) NF_STRERROR(ierr)
773            stop "getvarup"
774         endif
775!          write(*,*)'lecture dq ok',dq
776
777#ifdef NC_DOUBLE
778         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(23),hq)
779#else
780         ierr = NF_GET_VAR_REAL(nid,var3didin(23),hq)
781#endif
782         if(ierr/=NF_NOERR) then
783            write(*,*) NF_STRERROR(ierr)
784            stop "getvarup"
785         endif
786!          write(*,*)'lecture hq ok',hq
787
788#ifdef NC_DOUBLE
789         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(24),vq)
790#else
791         ierr = NF_GET_VAR_REAL(nid,var3didin(24),vq)
792#endif
793         if(ierr/=NF_NOERR) then
794            write(*,*) NF_STRERROR(ierr)
795            stop "getvarup"
796         endif
797!          write(*,*)'lecture vq ok',vq
798
799#ifdef NC_DOUBLE
800         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(25),dth)
801#else
802         ierr = NF_GET_VAR_REAL(nid,var3didin(25),dth)
803#endif
804         if(ierr/=NF_NOERR) then
805            write(*,*) NF_STRERROR(ierr)
806            stop "getvarup"
807         endif
808!          write(*,*)'lecture dth ok',dth
809
810#ifdef NC_DOUBLE
811         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(26),hth)
812#else
813         ierr = NF_GET_VAR_REAL(nid,var3didin(26),hth)
814#endif
815         if(ierr/=NF_NOERR) then
816            write(*,*) NF_STRERROR(ierr)
817            stop "getvarup"
818         endif
819!          write(*,*)'lecture hth ok',hth
820
821#ifdef NC_DOUBLE
822         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(27),vth)
823#else
824         ierr = NF_GET_VAR_REAL(nid,var3didin(27),vth)
825#endif
826         if(ierr/=NF_NOERR) then
827            write(*,*) NF_STRERROR(ierr)
828            stop "getvarup"
829         endif
830!          write(*,*)'lecture vth ok',vth
831
832#ifdef NC_DOUBLE
833         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(28),dr)
834#else
835         ierr = NF_GET_VAR_REAL(nid,var3didin(28),dr)
836#endif
837         if(ierr/=NF_NOERR) then
838            write(*,*) NF_STRERROR(ierr)
839            stop "getvarup"
840         endif
841!          write(*,*)'lecture dr ok',dr
842
843#ifdef NC_DOUBLE
844         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(29),hr)
845#else
846         ierr = NF_GET_VAR_REAL(nid,var3didin(29),hr)
847#endif
848         if(ierr/=NF_NOERR) then
849            write(*,*) NF_STRERROR(ierr)
850            stop "getvarup"
851         endif
852!          write(*,*)'lecture hr ok',hr
853
854#ifdef NC_DOUBLE
855         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(30),vr)
856#else
857         ierr = NF_GET_VAR_REAL(nid,var3didin(30),vr)
858#endif
859         if(ierr/=NF_NOERR) then
860            write(*,*) NF_STRERROR(ierr)
861            stop "getvarup"
862         endif
863!          write(*,*)'lecture vr ok',vr
864
865#ifdef NC_DOUBLE
866         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(31),dtrad)
867#else
868         ierr = NF_GET_VAR_REAL(nid,var3didin(31),dtrad)
869#endif
870         if(ierr/=NF_NOERR) then
871            write(*,*) NF_STRERROR(ierr)
872            stop "getvarup"
873         endif
874!          write(*,*)'lecture dtrad ok',dtrad
875
876#ifdef NC_DOUBLE
877         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(32),sens)
878#else
879         ierr = NF_GET_VAR_REAL(nid,var3didin(32),sens)
880#endif
881         if(ierr/=NF_NOERR) then
882            write(*,*) NF_STRERROR(ierr)
883            stop "getvarup"
884         endif
885!          write(*,*)'lecture sens ok',sens
886
887#ifdef NC_DOUBLE
888         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(33),flat)
889#else
890         ierr = NF_GET_VAR_REAL(nid,var3didin(33),flat)
891#endif
892         if(ierr/=NF_NOERR) then
893            write(*,*) NF_STRERROR(ierr)
894            stop "getvarup"
895         endif
896!          write(*,*)'lecture flat ok',flat
897
898#ifdef NC_DOUBLE
899         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(34),ts)
900#else
901         ierr = NF_GET_VAR_REAL(nid,var3didin(34),ts)
902#endif
903         if(ierr/=NF_NOERR) then
904            write(*,*) NF_STRERROR(ierr)
905            stop "getvarup"
906         endif
907!          write(*,*)'lecture ts ok',ts
908
909#ifdef NC_DOUBLE
910         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(35),ustar)
911#else
912         ierr = NF_GET_VAR_REAL(nid,var3didin(35),ustar)
913#endif
914         if(ierr/=NF_NOERR) then
915            write(*,*) NF_STRERROR(ierr)
916            stop "getvarup"
917         endif
918!         write(*,*)'lecture ustar ok',ustar
919
920#ifdef NC_DOUBLE
921         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(36),uw)
922#else
923         ierr = NF_GET_VAR_REAL(nid,var3didin(36),uw)
924#endif
925         if(ierr/=NF_NOERR) then
926            write(*,*) NF_STRERROR(ierr)
927            stop "getvarup"
928         endif
929!         write(*,*)'lecture uw ok',uw
930
931#ifdef NC_DOUBLE
932         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(37),vw)
933#else
934         ierr = NF_GET_VAR_REAL(nid,var3didin(37),vw)
935#endif
936         if(ierr/=NF_NOERR) then
937            write(*,*) NF_STRERROR(ierr)
938            stop "getvarup"
939         endif
940!         write(*,*)'lecture vw ok',vw
941
942#ifdef NC_DOUBLE
943         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(38),q1)
944#else
945         ierr = NF_GET_VAR_REAL(nid,var3didin(38),q1)
946#endif
947         if(ierr/=NF_NOERR) then
948            write(*,*) NF_STRERROR(ierr)
949            stop "getvarup"
950         endif
951!         write(*,*)'lecture q1 ok',q1
952
953#ifdef NC_DOUBLE
954         ierr = NF_GET_VAR_DOUBLE(nid,var3didin(39),q2)
955#else
956         ierr = NF_GET_VAR_REAL(nid,var3didin(39),q2)
957#endif
958         if(ierr/=NF_NOERR) then
959            write(*,*) NF_STRERROR(ierr)
960            stop "getvarup"
961         endif
962!         write(*,*)'lecture q2 ok',q2
963
964
965         return
966         end subroutine read_cas
967!======================================================================
968        SUBROUTINE interp_case_time(day,day1,annee_ref                &
969!    &         ,year_cas,day_cas,nt_cas,pdt_forc,nlev_cas      &
970     &         ,nt_cas,nlev_cas                                       &
971     &         ,ts_cas,plev_cas,t_cas,q_cas,u_cas,v_cas               &
972     &         ,ug_cas,vg_cas,vitw_cas,du_cas,hu_cas,vu_cas           &
973     &         ,dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dtrad_cas   &
974     &         ,dq_cas,hq_cas,vq_cas,lat_cas,sens_cas,ustar_cas       &
975     &         ,uw_cas,vw_cas,q1_cas,q2_cas                           &
976     &         ,ts_prof_cas,plev_prof_cas,t_prof_cas,q_prof_cas       &
977     &         ,u_prof_cas,v_prof_cas,ug_prof_cas,vg_prof_cas         &
978     &         ,vitw_prof_cas,du_prof_cas,hu_prof_cas,vu_prof_cas     &
979     &         ,dv_prof_cas,hv_prof_cas,vv_prof_cas,dt_prof_cas       &
980     &         ,ht_prof_cas,vt_prof_cas,dtrad_prof_cas,dq_prof_cas    &
981     &         ,hq_prof_cas,vq_prof_cas,lat_prof_cas,sens_prof_cas    &
982     &         ,ustar_prof_cas,uw_prof_cas,vw_prof_cas,q1_prof_cas,q2_prof_cas)
983         
984
985        implicit none
986
987!---------------------------------------------------------------------------------------
988! Time interpolation of a 2D field to the timestep corresponding to day
989!
990! day: current julian day (e.g. 717538.2)
991! day1: first day of the simulation
992! nt_cas: total nb of data in the forcing
993! pdt_cas: total time interval (in sec) between 2 forcing data
994!---------------------------------------------------------------------------------------
995
996#include "compar1d.h"
997#include "date_cas.h"
998
999! inputs:
1000        integer annee_ref
1001        integer nt_cas,nlev_cas
1002        real day, day1,day_cas
1003        real ts_cas(nt_cas)
1004        real plev_cas(nlev_cas,nt_cas)
1005        real t_cas(nlev_cas,nt_cas),q_cas(nlev_cas,nt_cas)
1006        real u_cas(nlev_cas,nt_cas),v_cas(nlev_cas,nt_cas)
1007        real ug_cas(nlev_cas,nt_cas),vg_cas(nlev_cas,nt_cas)
1008        real vitw_cas(nlev_cas,nt_cas)
1009        real du_cas(nlev_cas,nt_cas),hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas)
1010        real dv_cas(nlev_cas,nt_cas),hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas)
1011        real dt_cas(nlev_cas,nt_cas),ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas)
1012        real dtrad_cas(nlev_cas,nt_cas)
1013        real dq_cas(nlev_cas,nt_cas),hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas)
1014        real lat_cas(nt_cas)
1015        real sens_cas(nt_cas)
1016        real ustar_cas(nt_cas),uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas)
1017        real q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas)
1018
1019! outputs:
1020        real plev_prof_cas(nlev_cas)
1021        real t_prof_cas(nlev_cas),q_prof_cas(nlev_cas)
1022        real u_prof_cas(nlev_cas),v_prof_cas(nlev_cas)
1023        real ug_prof_cas(nlev_cas),vg_prof_cas(nlev_cas)
1024        real vitw_prof_cas(nlev_cas)
1025        real du_prof_cas(nlev_cas),hu_prof_cas(nlev_cas),vu_prof_cas(nlev_cas)
1026        real dv_prof_cas(nlev_cas),hv_prof_cas(nlev_cas),vv_prof_cas(nlev_cas)
1027        real dt_prof_cas(nlev_cas),ht_prof_cas(nlev_cas),vt_prof_cas(nlev_cas)
1028        real dtrad_prof_cas(nlev_cas)
1029        real dq_prof_cas(nlev_cas),hq_prof_cas(nlev_cas),vq_prof_cas(nlev_cas)
1030        real lat_prof_cas,sens_prof_cas,ts_prof_cas,ustar_prof_cas
1031        real uw_prof_cas(nlev_cas),vw_prof_cas(nlev_cas),q1_prof_cas(nlev_cas),q2_prof_cas(nlev_cas)
1032! local:
1033        integer it_cas1, it_cas2,k
1034        real timeit,time_cas1,time_cas2,frac
1035
1036
1037        print*,'Check time',day1,day_ju_ini_cas,day_deb+1,pdt_cas
1038
1039! On teste si la date du cas AMMA est correcte.
1040! C est pour memoire car en fait les fichiers .def
1041! sont censes etre corrects.
1042! A supprimer a terme (MPL 20150623)
1043!     if ((forcing_type.eq.10).and.(1.eq.0)) then
1044! Check that initial day of the simulation consistent with AMMA case:
1045!      if (annee_ref.ne.2006) then
1046!       print*,'Pour AMMA, annee_ref doit etre 2006'
1047!       print*,'Changer annee_ref dans run.def'
1048!       stop
1049!      endif
1050!      if (annee_ref.eq.2006 .and. day1.lt.day_cas) then
1051!       print*,'AMMA a debute le 10 juillet 2006',day1,day_cas
1052!       print*,'Changer dayref dans run.def'
1053!       stop
1054!      endif
1055!      if (annee_ref.eq.2006 .and. day1.gt.day_cas+1) then
1056!       print*,'AMMA a fini le 11 juillet'
1057!       print*,'Changer dayref ou nday dans run.def'
1058!       stop
1059!      endif
1060!      endif
1061
1062! Determine timestep relative to the 1st day:
1063!       timeit=(day-day1)*86400.
1064!       if (annee_ref.eq.1992) then
1065!        timeit=(day-day_cas)*86400.
1066!       else
1067!        timeit=(day+61.-1.)*86400. ! 61 days between Nov01 and Dec31 1992
1068!       endif
1069      timeit=(day-day_ju_ini_cas)*86400
1070      print *,'day=',day
1071      print *,'day_ju_ini_cas=',day_ju_ini_cas
1072      print *,'pdt_cas=',pdt_cas
1073      print *,'timeit=',timeit
1074      print *,'nt_cas=',nt_cas
1075
1076! Determine the closest observation times:
1077!       it_cas1=INT(timeit/pdt_cas)+1
1078!       it_cas2=it_cas1 + 1
1079!       time_cas1=(it_cas1-1)*pdt_cas
1080!       time_cas2=(it_cas2-1)*pdt_cas
1081
1082       it_cas1=INT(timeit/pdt_cas)+1
1083       IF (it_cas1 .EQ. nt_cas) THEN
1084       it_cas2=it_cas1
1085       ELSE
1086       it_cas2=it_cas1 + 1
1087       ENDIF
1088       time_cas1=(it_cas1-1)*pdt_cas
1089       time_cas2=(it_cas2-1)*pdt_cas
1090      print *,'it_cas1=',it_cas1
1091      print *,'it_cas2=',it_cas2
1092      print *,'time_cas1=',time_cas1
1093      print *,'time_cas2=',time_cas2
1094
1095       if (it_cas1 .gt. nt_cas) then
1096        write(*,*) 'PB-stop: day, day_ju_ini_cas,it_cas1, it_cas2, timeit: '            &
1097     &        ,day,day_ju_ini_cas,it_cas1,it_cas2,timeit
1098        stop
1099       endif
1100
1101! time interpolation:
1102       IF (it_cas1 .EQ. it_cas2) THEN
1103          frac=0.
1104       ELSE
1105          frac=(time_cas2-timeit)/(time_cas2-time_cas1)
1106          frac=max(frac,0.0)
1107       ENDIF
1108
1109       lat_prof_cas = lat_cas(it_cas2)                                       &
1110     &          -frac*(lat_cas(it_cas2)-lat_cas(it_cas1))
1111       sens_prof_cas = sens_cas(it_cas2)                                     &
1112     &          -frac*(sens_cas(it_cas2)-sens_cas(it_cas1))
1113       ts_prof_cas = ts_cas(it_cas2)                                         &
1114     &          -frac*(ts_cas(it_cas2)-ts_cas(it_cas1))
1115       ustar_prof_cas = ustar_cas(it_cas2)                                   &
1116     &          -frac*(ustar_cas(it_cas2)-ustar_cas(it_cas1))
1117
1118       do k=1,nlev_cas
1119        plev_prof_cas(k) = plev_cas(k,it_cas2)                               &
1120     &          -frac*(plev_cas(k,it_cas2)-plev_cas(k,it_cas1))
1121        t_prof_cas(k) = t_cas(k,it_cas2)                               &
1122     &          -frac*(t_cas(k,it_cas2)-t_cas(k,it_cas1))
1123        q_prof_cas(k) = q_cas(k,it_cas2)                               &
1124     &          -frac*(q_cas(k,it_cas2)-q_cas(k,it_cas1))
1125        u_prof_cas(k) = u_cas(k,it_cas2)                               &
1126     &          -frac*(u_cas(k,it_cas2)-u_cas(k,it_cas1))
1127        v_prof_cas(k) = v_cas(k,it_cas2)                               &
1128     &          -frac*(v_cas(k,it_cas2)-v_cas(k,it_cas1))
1129        ug_prof_cas(k) = ug_cas(k,it_cas2)                               &
1130     &          -frac*(ug_cas(k,it_cas2)-ug_cas(k,it_cas1))
1131        vg_prof_cas(k) = vg_cas(k,it_cas2)                               &
1132     &          -frac*(vg_cas(k,it_cas2)-vg_cas(k,it_cas1))
1133        vitw_prof_cas(k) = vitw_cas(k,it_cas2)                               &
1134     &          -frac*(vitw_cas(k,it_cas2)-vitw_cas(k,it_cas1))
1135        du_prof_cas(k) = du_cas(k,it_cas2)                                   &
1136     &          -frac*(du_cas(k,it_cas2)-du_cas(k,it_cas1))
1137        hu_prof_cas(k) = hu_cas(k,it_cas2)                                   &
1138     &          -frac*(hu_cas(k,it_cas2)-hu_cas(k,it_cas1))
1139        vu_prof_cas(k) = vu_cas(k,it_cas2)                                   &
1140     &          -frac*(vu_cas(k,it_cas2)-vu_cas(k,it_cas1))
1141        dv_prof_cas(k) = dv_cas(k,it_cas2)                                   &
1142     &          -frac*(dv_cas(k,it_cas2)-dv_cas(k,it_cas1))
1143        hv_prof_cas(k) = hv_cas(k,it_cas2)                                   &
1144     &          -frac*(hv_cas(k,it_cas2)-hv_cas(k,it_cas1))
1145        vv_prof_cas(k) = vv_cas(k,it_cas2)                                   &
1146     &          -frac*(vv_cas(k,it_cas2)-vv_cas(k,it_cas1))
1147        dt_prof_cas(k) = dt_cas(k,it_cas2)                                   &
1148     &          -frac*(dt_cas(k,it_cas2)-dt_cas(k,it_cas1))
1149        ht_prof_cas(k) = ht_cas(k,it_cas2)                                   &
1150     &          -frac*(ht_cas(k,it_cas2)-ht_cas(k,it_cas1))
1151        vt_prof_cas(k) = vt_cas(k,it_cas2)                                   &
1152     &          -frac*(vt_cas(k,it_cas2)-vt_cas(k,it_cas1))
1153        dtrad_prof_cas(k) = dtrad_cas(k,it_cas2)                                   &
1154     &          -frac*(dtrad_cas(k,it_cas2)-dtrad_cas(k,it_cas1))
1155        dq_prof_cas(k) = dq_cas(k,it_cas2)                                   &
1156     &          -frac*(dq_cas(k,it_cas2)-dq_cas(k,it_cas1))
1157        hq_prof_cas(k) = hq_cas(k,it_cas2)                                   &
1158     &          -frac*(hq_cas(k,it_cas2)-hq_cas(k,it_cas1))
1159        vq_prof_cas(k) = vq_cas(k,it_cas2)                                   &
1160     &          -frac*(vq_cas(k,it_cas2)-vq_cas(k,it_cas1))
1161       uw_prof_cas(k) = uw_cas(k,it_cas2)                                   &
1162     &          -frac*(uw_cas(k,it_cas2)-uw_cas(k,it_cas1))
1163       vw_prof_cas(k) = vw_cas(k,it_cas2)                                   &
1164     &          -frac*(vw_cas(k,it_cas2)-vw_cas(k,it_cas1))
1165       q1_prof_cas(k) = q1_cas(k,it_cas2)                                   &
1166     &          -frac*(q1_cas(k,it_cas2)-q1_cas(k,it_cas1))
1167       q2_prof_cas(k) = q2_cas(k,it_cas2)                                   &
1168     &          -frac*(q2_cas(k,it_cas2)-q2_cas(k,it_cas1))
1169        enddo
1170
1171        return
1172        END
1173
1174!**********************************************************************************************
Note: See TracBrowser for help on using the repository browser.