source: LMDZ6/branches/Amaury_dev/libf/phylmd/dyn1d/mod_1D_cases_read.F90

Last change on this file was 5295, checked in by abarral, 2 months ago

lint
F90 -> f90
finish removing svn !Id: headers, which were inconsistent

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