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

Last change on this file since 5260 was 5249, checked in by abarral, 5 weeks ago

Replace uses of cpp key NC_DOUBLE

  • Property svn:keywords set to Id
File size: 34.7 KB
Line 
1!
2! $Id: mod_1D_cases_read.F90 5249 2024-10-22 09:35:08Z abarral $
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      USE netcdf, ONLY: nf90_get_var
269      implicit none
270      INCLUDE "netcdf.inc"
271
272      integer ntime,nlevel
273
274      real zz(nlevel,ntime)
275      real pp(nlevel,ntime)
276      real temp(nlevel,ntime),qv(nlevel,ntime),rh(nlevel,ntime)
277      real theta(nlevel,ntime),rv(nlevel,ntime)
278      real u(nlevel,ntime)
279      real v(nlevel,ntime)
280      real ug(nlevel,ntime)
281      real vg(nlevel,ntime)
282      real w(nlevel,ntime)
283      real du(nlevel,ntime),hu(nlevel,ntime),vu(nlevel,ntime)
284      real dv(nlevel,ntime),hv(nlevel,ntime),vv(nlevel,ntime)
285      real dt(nlevel,ntime),ht(nlevel,ntime),vt(nlevel,ntime)
286      real dtrad(nlevel,ntime)
287      real dq(nlevel,ntime),hq(nlevel,ntime),vq(nlevel,ntime)
288      real dth(nlevel,ntime),hth(nlevel,ntime),vth(nlevel,ntime)
289      real dr(nlevel,ntime),hr(nlevel,ntime),vr(nlevel,ntime)
290      real flat(ntime),sens(ntime),ts(ntime),ustar(ntime)
291      real uw(nlevel,ntime),vw(nlevel,ntime),q1(nlevel,ntime),q2(nlevel,ntime)
292
293
294      integer nid, ierr,rid
295      integer nbvar3d
296      parameter(nbvar3d=39)
297      integer var3didin(nbvar3d)
298
299       ierr=NF_INQ_VARID(nid,"zz",var3didin(1))
300         if(ierr/=NF_NOERR) then
301           write(*,*) NF_STRERROR(ierr)
302           stop 'lev'
303         endif
304     
305      ierr=NF_INQ_VARID(nid,"pp",var3didin(2))
306         if(ierr/=NF_NOERR) then
307           write(*,*) NF_STRERROR(ierr)
308           stop 'plev'
309         endif
310
311
312      ierr=NF_INQ_VARID(nid,"temp",var3didin(3))
313         if(ierr/=NF_NOERR) then
314           write(*,*) NF_STRERROR(ierr)
315           stop 'temp'
316         endif
317
318      ierr=NF_INQ_VARID(nid,"qv",var3didin(4))
319         if(ierr/=NF_NOERR) then
320           write(*,*) NF_STRERROR(ierr)
321           stop 'qv'
322         endif
323
324      ierr=NF_INQ_VARID(nid,"rh",var3didin(5))
325         if(ierr/=NF_NOERR) then
326           write(*,*) NF_STRERROR(ierr)
327           stop 'rh'
328         endif
329
330      ierr=NF_INQ_VARID(nid,"theta",var3didin(6))
331         if(ierr/=NF_NOERR) then
332           write(*,*) NF_STRERROR(ierr)
333           stop 'theta'
334         endif
335
336      ierr=NF_INQ_VARID(nid,"rv",var3didin(7))
337         if(ierr/=NF_NOERR) then
338           write(*,*) NF_STRERROR(ierr)
339           stop 'rv'
340         endif
341
342
343      ierr=NF_INQ_VARID(nid,"u",var3didin(8))
344         if(ierr/=NF_NOERR) then
345           write(*,*) NF_STRERROR(ierr)
346           stop 'u'
347         endif
348
349      ierr=NF_INQ_VARID(nid,"v",var3didin(9))
350         if(ierr/=NF_NOERR) then
351           write(*,*) NF_STRERROR(ierr)
352           stop 'v'
353         endif
354
355       ierr=NF_INQ_VARID(nid,"ug",var3didin(10))
356         if(ierr/=NF_NOERR) then
357           write(*,*) NF_STRERROR(ierr)
358           stop 'ug'
359         endif
360
361      ierr=NF_INQ_VARID(nid,"vg",var3didin(11))
362         if(ierr/=NF_NOERR) then
363           write(*,*) NF_STRERROR(ierr)
364           stop 'vg'
365         endif
366
367      ierr=NF_INQ_VARID(nid,"w",var3didin(12))
368         if(ierr/=NF_NOERR) then
369           write(*,*) NF_STRERROR(ierr)
370           stop 'w'
371         endif
372
373      ierr=NF_INQ_VARID(nid,"advu",var3didin(13))
374         if(ierr/=NF_NOERR) then
375           write(*,*) NF_STRERROR(ierr)
376           stop 'advu'
377         endif
378
379      ierr=NF_INQ_VARID(nid,"hu",var3didin(14))
380         if(ierr/=NF_NOERR) then
381           write(*,*) NF_STRERROR(ierr)
382           stop 'hu'
383         endif
384
385       ierr=NF_INQ_VARID(nid,"vu",var3didin(15))
386         if(ierr/=NF_NOERR) then
387           write(*,*) NF_STRERROR(ierr)
388           stop 'vu'
389         endif
390
391       ierr=NF_INQ_VARID(nid,"advv",var3didin(16))
392         if(ierr/=NF_NOERR) then
393           write(*,*) NF_STRERROR(ierr)
394           stop 'advv'
395         endif
396
397      ierr=NF_INQ_VARID(nid,"hv",var3didin(17))
398         if(ierr/=NF_NOERR) then
399           write(*,*) NF_STRERROR(ierr)
400           stop 'hv'
401         endif
402
403       ierr=NF_INQ_VARID(nid,"vv",var3didin(18))
404         if(ierr/=NF_NOERR) then
405           write(*,*) NF_STRERROR(ierr)
406           stop 'vv'
407         endif
408
409      ierr=NF_INQ_VARID(nid,"advT",var3didin(19))
410         if(ierr/=NF_NOERR) then
411           write(*,*) NF_STRERROR(ierr)
412           stop 'advT'
413         endif
414
415      ierr=NF_INQ_VARID(nid,"hT",var3didin(20))
416         if(ierr/=NF_NOERR) then
417           write(*,*) NF_STRERROR(ierr)
418           stop 'hT'
419         endif
420
421      ierr=NF_INQ_VARID(nid,"vT",var3didin(21))
422         if(ierr/=NF_NOERR) then
423           write(*,*) NF_STRERROR(ierr)
424           stop 'vT'
425         endif
426
427      ierr=NF_INQ_VARID(nid,"advq",var3didin(22))
428         if(ierr/=NF_NOERR) then
429           write(*,*) NF_STRERROR(ierr)
430           stop 'advq'
431         endif
432     
433      ierr=NF_INQ_VARID(nid,"hq",var3didin(23))
434         if(ierr/=NF_NOERR) then
435           write(*,*) NF_STRERROR(ierr)
436           stop 'hq'
437         endif
438
439      ierr=NF_INQ_VARID(nid,"vq",var3didin(24))
440         if(ierr/=NF_NOERR) then
441           write(*,*) NF_STRERROR(ierr)
442           stop 'vq'
443         endif
444
445      ierr=NF_INQ_VARID(nid,"advth",var3didin(25))
446         if(ierr/=NF_NOERR) then
447           write(*,*) NF_STRERROR(ierr)
448           stop 'advth'
449         endif
450
451      ierr=NF_INQ_VARID(nid,"hth",var3didin(26))
452         if(ierr/=NF_NOERR) then
453           write(*,*) NF_STRERROR(ierr)
454           stop 'hth'
455         endif
456
457      ierr=NF_INQ_VARID(nid,"vth",var3didin(27))
458         if(ierr/=NF_NOERR) then
459           write(*,*) NF_STRERROR(ierr)
460           stop 'vth'
461         endif
462
463      ierr=NF_INQ_VARID(nid,"advr",var3didin(28))
464         if(ierr/=NF_NOERR) then
465           write(*,*) NF_STRERROR(ierr)
466           stop 'advr'
467         endif
468     
469      ierr=NF_INQ_VARID(nid,"hr",var3didin(29))
470         if(ierr/=NF_NOERR) then
471           write(*,*) NF_STRERROR(ierr)
472           stop 'hr'
473         endif
474
475      ierr=NF_INQ_VARID(nid,"vr",var3didin(30))
476         if(ierr/=NF_NOERR) then
477           write(*,*) NF_STRERROR(ierr)
478           stop 'vr'
479         endif
480
481      ierr=NF_INQ_VARID(nid,"radT",var3didin(31))
482         if(ierr/=NF_NOERR) then
483           write(*,*) NF_STRERROR(ierr)
484           stop 'radT'
485         endif
486
487      ierr=NF_INQ_VARID(nid,"sens",var3didin(32))
488         if(ierr/=NF_NOERR) then
489           write(*,*) NF_STRERROR(ierr)
490           stop 'sens'
491         endif
492
493      ierr=NF_INQ_VARID(nid,"flat",var3didin(33))
494         if(ierr/=NF_NOERR) then
495           write(*,*) NF_STRERROR(ierr)
496           stop 'flat'
497         endif
498
499      ierr=NF_INQ_VARID(nid,"ts",var3didin(34))
500         if(ierr/=NF_NOERR) then
501           write(*,*) NF_STRERROR(ierr)
502           stop 'ts'
503         endif
504
505      ierr=NF_INQ_VARID(nid,"ustar",var3didin(35))
506         if(ierr/=NF_NOERR) then
507           write(*,*) NF_STRERROR(ierr)
508           stop 'ustar'
509         endif
510
511      ierr=NF_INQ_VARID(nid,"uw",var3didin(36))
512         if(ierr/=NF_NOERR) then
513           write(*,*) NF_STRERROR(ierr)
514           stop 'uw'
515         endif
516
517      ierr=NF_INQ_VARID(nid,"vw",var3didin(37))
518         if(ierr/=NF_NOERR) then
519           write(*,*) NF_STRERROR(ierr)
520           stop 'vw'
521         endif
522
523      ierr=NF_INQ_VARID(nid,"q1",var3didin(38))
524         if(ierr/=NF_NOERR) then
525           write(*,*) NF_STRERROR(ierr)
526           stop 'q1'
527         endif
528
529      ierr=NF_INQ_VARID(nid,"q2",var3didin(39))
530         if(ierr/=NF_NOERR) then
531           write(*,*) NF_STRERROR(ierr)
532           stop 'q2'
533         endif
534 
535         ierr = nf90_get_var(nid, var3didin(1), zz)
536         if(ierr/=NF_NOERR) then
537            write(*,*) NF_STRERROR(ierr)
538            stop "getvarup"
539         endif
540!          write(*,*)'lecture z ok',zz
541
542         ierr = nf90_get_var(nid, var3didin(2), pp)
543         if(ierr/=NF_NOERR) then
544            write(*,*) NF_STRERROR(ierr)
545            stop "getvarup"
546         endif
547!          write(*,*)'lecture pp ok',pp
548
549
550         ierr = nf90_get_var(nid, var3didin(3), temp)
551         if(ierr/=NF_NOERR) then
552            write(*,*) NF_STRERROR(ierr)
553            stop "getvarup"
554         endif
555!          write(*,*)'lecture T ok',temp
556
557         ierr = nf90_get_var(nid, var3didin(4), qv)
558         if(ierr/=NF_NOERR) then
559            write(*,*) NF_STRERROR(ierr)
560            stop "getvarup"
561         endif
562!          write(*,*)'lecture qv ok',qv
563 
564         ierr = nf90_get_var(nid, var3didin(5), rh)
565         if(ierr/=NF_NOERR) then
566            write(*,*) NF_STRERROR(ierr)
567            stop "getvarup"
568         endif
569!          write(*,*)'lecture rh ok',rh
570
571         ierr = nf90_get_var(nid, var3didin(6), theta)
572         if(ierr/=NF_NOERR) then
573            write(*,*) NF_STRERROR(ierr)
574            stop "getvarup"
575         endif
576!          write(*,*)'lecture theta ok',theta
577
578         ierr = nf90_get_var(nid, var3didin(7), rv)
579         if(ierr/=NF_NOERR) then
580            write(*,*) NF_STRERROR(ierr)
581            stop "getvarup"
582         endif
583!          write(*,*)'lecture rv ok',rv
584
585         ierr = nf90_get_var(nid, var3didin(8), u)
586         if(ierr/=NF_NOERR) then
587            write(*,*) NF_STRERROR(ierr)
588            stop "getvarup"
589         endif
590!          write(*,*)'lecture u ok',u
591
592         ierr = nf90_get_var(nid, var3didin(9), v)
593         if(ierr/=NF_NOERR) then
594            write(*,*) NF_STRERROR(ierr)
595            stop "getvarup"
596         endif
597!          write(*,*)'lecture v ok',v
598
599         ierr = nf90_get_var(nid, var3didin(10), ug)
600         if(ierr/=NF_NOERR) then
601            write(*,*) NF_STRERROR(ierr)
602            stop "getvarup"
603         endif
604!          write(*,*)'lecture ug ok',ug
605
606         ierr = nf90_get_var(nid, var3didin(11), vg)
607         if(ierr/=NF_NOERR) then
608            write(*,*) NF_STRERROR(ierr)
609            stop "getvarup"
610         endif
611!          write(*,*)'lecture vg ok',vg
612
613         ierr = nf90_get_var(nid, var3didin(12), w)
614         if(ierr/=NF_NOERR) then
615            write(*,*) NF_STRERROR(ierr)
616            stop "getvarup"
617         endif
618!          write(*,*)'lecture w ok',w
619
620         ierr = nf90_get_var(nid, var3didin(13), du)
621         if(ierr/=NF_NOERR) then
622            write(*,*) NF_STRERROR(ierr)
623            stop "getvarup"
624         endif
625!          write(*,*)'lecture du ok',du
626
627         ierr = nf90_get_var(nid, var3didin(14), hu)
628         if(ierr/=NF_NOERR) then
629            write(*,*) NF_STRERROR(ierr)
630            stop "getvarup"
631         endif
632!          write(*,*)'lecture hu ok',hu
633
634         ierr = nf90_get_var(nid, var3didin(15), vu)
635         if(ierr/=NF_NOERR) then
636            write(*,*) NF_STRERROR(ierr)
637            stop "getvarup"
638         endif
639!          write(*,*)'lecture vu ok',vu
640
641         ierr = nf90_get_var(nid, var3didin(16), dv)
642         if(ierr/=NF_NOERR) then
643            write(*,*) NF_STRERROR(ierr)
644            stop "getvarup"
645         endif
646!          write(*,*)'lecture dv ok',dv
647
648         ierr = nf90_get_var(nid, var3didin(17), hv)
649         if(ierr/=NF_NOERR) then
650            write(*,*) NF_STRERROR(ierr)
651            stop "getvarup"
652         endif
653!          write(*,*)'lecture hv ok',hv
654
655         ierr = nf90_get_var(nid, var3didin(18), vv)
656         if(ierr/=NF_NOERR) then
657            write(*,*) NF_STRERROR(ierr)
658            stop "getvarup"
659         endif
660!          write(*,*)'lecture vv ok',vv
661
662         ierr = nf90_get_var(nid, var3didin(19), dt)
663         if(ierr/=NF_NOERR) then
664            write(*,*) NF_STRERROR(ierr)
665            stop "getvarup"
666         endif
667!          write(*,*)'lecture dt ok',dt
668
669         ierr = nf90_get_var(nid, var3didin(20), ht)
670         if(ierr/=NF_NOERR) then
671            write(*,*) NF_STRERROR(ierr)
672            stop "getvarup"
673         endif
674!          write(*,*)'lecture ht ok',ht
675
676         ierr = nf90_get_var(nid, var3didin(21), vt)
677         if(ierr/=NF_NOERR) then
678            write(*,*) NF_STRERROR(ierr)
679            stop "getvarup"
680         endif
681!          write(*,*)'lecture vt ok',vt
682
683         ierr = nf90_get_var(nid, var3didin(22), dq)
684         if(ierr/=NF_NOERR) then
685            write(*,*) NF_STRERROR(ierr)
686            stop "getvarup"
687         endif
688!          write(*,*)'lecture dq ok',dq
689
690         ierr = nf90_get_var(nid, var3didin(23), hq)
691         if(ierr/=NF_NOERR) then
692            write(*,*) NF_STRERROR(ierr)
693            stop "getvarup"
694         endif
695!          write(*,*)'lecture hq ok',hq
696
697         ierr = nf90_get_var(nid, var3didin(24), vq)
698         if(ierr/=NF_NOERR) then
699            write(*,*) NF_STRERROR(ierr)
700            stop "getvarup"
701         endif
702!          write(*,*)'lecture vq ok',vq
703
704         ierr = nf90_get_var(nid, var3didin(25), dth)
705         if(ierr/=NF_NOERR) then
706            write(*,*) NF_STRERROR(ierr)
707            stop "getvarup"
708         endif
709!          write(*,*)'lecture dth ok',dth
710
711         ierr = nf90_get_var(nid, var3didin(26), hth)
712         if(ierr/=NF_NOERR) then
713            write(*,*) NF_STRERROR(ierr)
714            stop "getvarup"
715         endif
716!          write(*,*)'lecture hth ok',hth
717
718         ierr = nf90_get_var(nid, var3didin(27), vth)
719         if(ierr/=NF_NOERR) then
720            write(*,*) NF_STRERROR(ierr)
721            stop "getvarup"
722         endif
723!          write(*,*)'lecture vth ok',vth
724
725         ierr = nf90_get_var(nid, var3didin(28), dr)
726         if(ierr/=NF_NOERR) then
727            write(*,*) NF_STRERROR(ierr)
728            stop "getvarup"
729         endif
730!          write(*,*)'lecture dr ok',dr
731
732         ierr = nf90_get_var(nid, var3didin(29), hr)
733         if(ierr/=NF_NOERR) then
734            write(*,*) NF_STRERROR(ierr)
735            stop "getvarup"
736         endif
737!          write(*,*)'lecture hr ok',hr
738
739         ierr = nf90_get_var(nid, var3didin(30), vr)
740         if(ierr/=NF_NOERR) then
741            write(*,*) NF_STRERROR(ierr)
742            stop "getvarup"
743         endif
744!          write(*,*)'lecture vr ok',vr
745
746         ierr = nf90_get_var(nid, var3didin(31), dtrad)
747         if(ierr/=NF_NOERR) then
748            write(*,*) NF_STRERROR(ierr)
749            stop "getvarup"
750         endif
751!          write(*,*)'lecture dtrad ok',dtrad
752
753         ierr = nf90_get_var(nid, var3didin(32), sens)
754         if(ierr/=NF_NOERR) then
755            write(*,*) NF_STRERROR(ierr)
756            stop "getvarup"
757         endif
758!          write(*,*)'lecture sens ok',sens
759
760         ierr = nf90_get_var(nid, var3didin(33), flat)
761         if(ierr/=NF_NOERR) then
762            write(*,*) NF_STRERROR(ierr)
763            stop "getvarup"
764         endif
765!          write(*,*)'lecture flat ok',flat
766
767         ierr = nf90_get_var(nid, var3didin(34), ts)
768         if(ierr/=NF_NOERR) then
769            write(*,*) NF_STRERROR(ierr)
770            stop "getvarup"
771         endif
772!          write(*,*)'lecture ts ok',ts
773
774         ierr = nf90_get_var(nid, var3didin(35), ustar)
775         if(ierr/=NF_NOERR) then
776            write(*,*) NF_STRERROR(ierr)
777            stop "getvarup"
778         endif
779!         write(*,*)'lecture ustar ok',ustar
780
781         ierr = nf90_get_var(nid, var3didin(36), uw)
782         if(ierr/=NF_NOERR) then
783            write(*,*) NF_STRERROR(ierr)
784            stop "getvarup"
785         endif
786!         write(*,*)'lecture uw ok',uw
787
788         ierr = nf90_get_var(nid, var3didin(37), vw)
789         if(ierr/=NF_NOERR) then
790            write(*,*) NF_STRERROR(ierr)
791            stop "getvarup"
792         endif
793!         write(*,*)'lecture vw ok',vw
794
795         ierr = nf90_get_var(nid, var3didin(38), q1)
796         if(ierr/=NF_NOERR) then
797            write(*,*) NF_STRERROR(ierr)
798            stop "getvarup"
799         endif
800!         write(*,*)'lecture q1 ok',q1
801
802         ierr = nf90_get_var(nid, var3didin(39), q2)
803         if(ierr/=NF_NOERR) then
804            write(*,*) NF_STRERROR(ierr)
805            stop "getvarup"
806         endif
807!         write(*,*)'lecture q2 ok',q2
808
809
810         return
811         end subroutine read_cas
812!======================================================================
813        SUBROUTINE interp_case_time(day,day1,annee_ref                &
814!    &         ,year_cas,day_cas,nt_cas,pdt_forc,nlev_cas      &
815     &         ,nt_cas,nlev_cas                                       &
816     &         ,ts_cas,plev_cas,t_cas,q_cas,u_cas,v_cas               &
817     &         ,ug_cas,vg_cas,vitw_cas,du_cas,hu_cas,vu_cas           &
818     &         ,dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dtrad_cas   &
819     &         ,dq_cas,hq_cas,vq_cas,lat_cas,sens_cas,ustar_cas       &
820     &         ,uw_cas,vw_cas,q1_cas,q2_cas                           &
821     &         ,ts_prof_cas,plev_prof_cas,t_prof_cas,q_prof_cas       &
822     &         ,u_prof_cas,v_prof_cas,ug_prof_cas,vg_prof_cas         &
823     &         ,vitw_prof_cas,du_prof_cas,hu_prof_cas,vu_prof_cas     &
824     &         ,dv_prof_cas,hv_prof_cas,vv_prof_cas,dt_prof_cas       &
825     &         ,ht_prof_cas,vt_prof_cas,dtrad_prof_cas,dq_prof_cas    &
826     &         ,hq_prof_cas,vq_prof_cas,lat_prof_cas,sens_prof_cas    &
827     &         ,ustar_prof_cas,uw_prof_cas,vw_prof_cas,q1_prof_cas,q2_prof_cas)
828         
829
830        implicit none
831
832!---------------------------------------------------------------------------------------
833! Time interpolation of a 2D field to the timestep corresponding to day
834!
835! day: current julian day (e.g. 717538.2)
836! day1: first day of the simulation
837! nt_cas: total nb of data in the forcing
838! pdt_cas: total time interval (in sec) between 2 forcing data
839!---------------------------------------------------------------------------------------
840
841        INCLUDE "compar1d.h"
842        INCLUDE "date_cas.h"
843
844! inputs:
845        integer annee_ref
846        integer nt_cas,nlev_cas
847        real day, day1,day_cas
848        real ts_cas(nt_cas)
849        real plev_cas(nlev_cas,nt_cas)
850        real t_cas(nlev_cas,nt_cas),q_cas(nlev_cas,nt_cas)
851        real u_cas(nlev_cas,nt_cas),v_cas(nlev_cas,nt_cas)
852        real ug_cas(nlev_cas,nt_cas),vg_cas(nlev_cas,nt_cas)
853        real vitw_cas(nlev_cas,nt_cas)
854        real du_cas(nlev_cas,nt_cas),hu_cas(nlev_cas,nt_cas),vu_cas(nlev_cas,nt_cas)
855        real dv_cas(nlev_cas,nt_cas),hv_cas(nlev_cas,nt_cas),vv_cas(nlev_cas,nt_cas)
856        real dt_cas(nlev_cas,nt_cas),ht_cas(nlev_cas,nt_cas),vt_cas(nlev_cas,nt_cas)
857        real dtrad_cas(nlev_cas,nt_cas)
858        real dq_cas(nlev_cas,nt_cas),hq_cas(nlev_cas,nt_cas),vq_cas(nlev_cas,nt_cas)
859        real lat_cas(nt_cas)
860        real sens_cas(nt_cas)
861        real ustar_cas(nt_cas),uw_cas(nlev_cas,nt_cas),vw_cas(nlev_cas,nt_cas)
862        real q1_cas(nlev_cas,nt_cas),q2_cas(nlev_cas,nt_cas)
863
864! outputs:
865        real plev_prof_cas(nlev_cas)
866        real t_prof_cas(nlev_cas),q_prof_cas(nlev_cas)
867        real u_prof_cas(nlev_cas),v_prof_cas(nlev_cas)
868        real ug_prof_cas(nlev_cas),vg_prof_cas(nlev_cas)
869        real vitw_prof_cas(nlev_cas)
870        real du_prof_cas(nlev_cas),hu_prof_cas(nlev_cas),vu_prof_cas(nlev_cas)
871        real dv_prof_cas(nlev_cas),hv_prof_cas(nlev_cas),vv_prof_cas(nlev_cas)
872        real dt_prof_cas(nlev_cas),ht_prof_cas(nlev_cas),vt_prof_cas(nlev_cas)
873        real dtrad_prof_cas(nlev_cas)
874        real dq_prof_cas(nlev_cas),hq_prof_cas(nlev_cas),vq_prof_cas(nlev_cas)
875        real lat_prof_cas,sens_prof_cas,ts_prof_cas,ustar_prof_cas
876        real uw_prof_cas(nlev_cas),vw_prof_cas(nlev_cas),q1_prof_cas(nlev_cas),q2_prof_cas(nlev_cas)
877! local:
878        integer it_cas1, it_cas2,k
879        real timeit,time_cas1,time_cas2,frac
880
881
882        print*,'Check time',day1,day_ju_ini_cas,day_deb+1,pdt_cas
883
884! On teste si la date du cas AMMA est correcte.
885! C est pour memoire car en fait les fichiers .def
886! sont censes etre corrects.
887! A supprimer a terme (MPL 20150623)
888!     if ((forcing_type.eq.10).and.(1.eq.0)) then
889! Check that initial day of the simulation consistent with AMMA case:
890!      if (annee_ref.ne.2006) then
891!       print*,'Pour AMMA, annee_ref doit etre 2006'
892!       print*,'Changer annee_ref dans run.def'
893!       stop
894!      endif
895!      if (annee_ref.eq.2006 .and. day1.lt.day_cas) then
896!       print*,'AMMA a debute le 10 juillet 2006',day1,day_cas
897!       print*,'Changer dayref dans run.def'
898!       stop
899!      endif
900!      if (annee_ref.eq.2006 .and. day1.gt.day_cas+1) then
901!       print*,'AMMA a fini le 11 juillet'
902!       print*,'Changer dayref ou nday dans run.def'
903!       stop
904!      endif
905!      endif
906
907! Determine timestep relative to the 1st day:
908!       timeit=(day-day1)*86400.
909!       if (annee_ref.eq.1992) then
910!        timeit=(day-day_cas)*86400.
911!       else
912!        timeit=(day+61.-1.)*86400. ! 61 days between Nov01 and Dec31 1992
913!       endif
914      timeit=(day-day_ju_ini_cas)*86400
915      print *,'day=',day
916      print *,'day_ju_ini_cas=',day_ju_ini_cas
917      print *,'pdt_cas=',pdt_cas
918      print *,'timeit=',timeit
919      print *,'nt_cas=',nt_cas
920
921! Determine the closest observation times:
922!       it_cas1=INT(timeit/pdt_cas)+1
923!       it_cas2=it_cas1 + 1
924!       time_cas1=(it_cas1-1)*pdt_cas
925!       time_cas2=(it_cas2-1)*pdt_cas
926
927       it_cas1=INT(timeit/pdt_cas)+1
928       IF (it_cas1 .EQ. nt_cas) THEN
929       it_cas2=it_cas1
930       ELSE
931       it_cas2=it_cas1 + 1
932       ENDIF
933       time_cas1=(it_cas1-1)*pdt_cas
934       time_cas2=(it_cas2-1)*pdt_cas
935      print *,'it_cas1=',it_cas1
936      print *,'it_cas2=',it_cas2
937      print *,'time_cas1=',time_cas1
938      print *,'time_cas2=',time_cas2
939
940       if (it_cas1 .gt. nt_cas) then
941        write(*,*) 'PB-stop: day, day_ju_ini_cas,it_cas1, it_cas2, timeit: '            &
942     &        ,day,day_ju_ini_cas,it_cas1,it_cas2,timeit
943        stop
944       endif
945
946! time interpolation:
947       IF (it_cas1 .EQ. it_cas2) THEN
948          frac=0.
949       ELSE
950          frac=(time_cas2-timeit)/(time_cas2-time_cas1)
951          frac=max(frac,0.0)
952       ENDIF
953
954       lat_prof_cas = lat_cas(it_cas2)                                       &
955     &          -frac*(lat_cas(it_cas2)-lat_cas(it_cas1))
956       sens_prof_cas = sens_cas(it_cas2)                                     &
957     &          -frac*(sens_cas(it_cas2)-sens_cas(it_cas1))
958       ts_prof_cas = ts_cas(it_cas2)                                         &
959     &          -frac*(ts_cas(it_cas2)-ts_cas(it_cas1))
960       ustar_prof_cas = ustar_cas(it_cas2)                                   &
961     &          -frac*(ustar_cas(it_cas2)-ustar_cas(it_cas1))
962
963       do k=1,nlev_cas
964        plev_prof_cas(k) = plev_cas(k,it_cas2)                               &
965     &          -frac*(plev_cas(k,it_cas2)-plev_cas(k,it_cas1))
966        t_prof_cas(k) = t_cas(k,it_cas2)                               &
967     &          -frac*(t_cas(k,it_cas2)-t_cas(k,it_cas1))
968        q_prof_cas(k) = q_cas(k,it_cas2)                               &
969     &          -frac*(q_cas(k,it_cas2)-q_cas(k,it_cas1))
970        u_prof_cas(k) = u_cas(k,it_cas2)                               &
971     &          -frac*(u_cas(k,it_cas2)-u_cas(k,it_cas1))
972        v_prof_cas(k) = v_cas(k,it_cas2)                               &
973     &          -frac*(v_cas(k,it_cas2)-v_cas(k,it_cas1))
974        ug_prof_cas(k) = ug_cas(k,it_cas2)                               &
975     &          -frac*(ug_cas(k,it_cas2)-ug_cas(k,it_cas1))
976        vg_prof_cas(k) = vg_cas(k,it_cas2)                               &
977     &          -frac*(vg_cas(k,it_cas2)-vg_cas(k,it_cas1))
978        vitw_prof_cas(k) = vitw_cas(k,it_cas2)                               &
979     &          -frac*(vitw_cas(k,it_cas2)-vitw_cas(k,it_cas1))
980        du_prof_cas(k) = du_cas(k,it_cas2)                                   &
981     &          -frac*(du_cas(k,it_cas2)-du_cas(k,it_cas1))
982        hu_prof_cas(k) = hu_cas(k,it_cas2)                                   &
983     &          -frac*(hu_cas(k,it_cas2)-hu_cas(k,it_cas1))
984        vu_prof_cas(k) = vu_cas(k,it_cas2)                                   &
985     &          -frac*(vu_cas(k,it_cas2)-vu_cas(k,it_cas1))
986        dv_prof_cas(k) = dv_cas(k,it_cas2)                                   &
987     &          -frac*(dv_cas(k,it_cas2)-dv_cas(k,it_cas1))
988        hv_prof_cas(k) = hv_cas(k,it_cas2)                                   &
989     &          -frac*(hv_cas(k,it_cas2)-hv_cas(k,it_cas1))
990        vv_prof_cas(k) = vv_cas(k,it_cas2)                                   &
991     &          -frac*(vv_cas(k,it_cas2)-vv_cas(k,it_cas1))
992        dt_prof_cas(k) = dt_cas(k,it_cas2)                                   &
993     &          -frac*(dt_cas(k,it_cas2)-dt_cas(k,it_cas1))
994        ht_prof_cas(k) = ht_cas(k,it_cas2)                                   &
995     &          -frac*(ht_cas(k,it_cas2)-ht_cas(k,it_cas1))
996        vt_prof_cas(k) = vt_cas(k,it_cas2)                                   &
997     &          -frac*(vt_cas(k,it_cas2)-vt_cas(k,it_cas1))
998        dtrad_prof_cas(k) = dtrad_cas(k,it_cas2)                                   &
999     &          -frac*(dtrad_cas(k,it_cas2)-dtrad_cas(k,it_cas1))
1000        dq_prof_cas(k) = dq_cas(k,it_cas2)                                   &
1001     &          -frac*(dq_cas(k,it_cas2)-dq_cas(k,it_cas1))
1002        hq_prof_cas(k) = hq_cas(k,it_cas2)                                   &
1003     &          -frac*(hq_cas(k,it_cas2)-hq_cas(k,it_cas1))
1004        vq_prof_cas(k) = vq_cas(k,it_cas2)                                   &
1005     &          -frac*(vq_cas(k,it_cas2)-vq_cas(k,it_cas1))
1006       uw_prof_cas(k) = uw_cas(k,it_cas2)                                   &
1007     &          -frac*(uw_cas(k,it_cas2)-uw_cas(k,it_cas1))
1008       vw_prof_cas(k) = vw_cas(k,it_cas2)                                   &
1009     &          -frac*(vw_cas(k,it_cas2)-vw_cas(k,it_cas1))
1010       q1_prof_cas(k) = q1_cas(k,it_cas2)                                   &
1011     &          -frac*(q1_cas(k,it_cas2)-q1_cas(k,it_cas1))
1012       q2_prof_cas(k) = q2_cas(k,it_cas2)                                   &
1013     &          -frac*(q2_cas(k,it_cas2)-q2_cas(k,it_cas1))
1014        enddo
1015
1016        return
1017        END
1018
1019!**********************************************************************************************
Note: See TracBrowser for help on using the repository browser.