source: LMDZ5/trunk/libf/phylmd/1D_read_forc_cases.h @ 2210

Last change on this file since 2210 was 2191, checked in by fhourdin, 10 years ago

Ajout du cas 1D CINDY-DYNAMO, utilisant le nouveau format standard, amené à être étendu aux autres cas.
Addition of the CINDY-DYNAMO 1D case, using the new standard format for 1D cases, that will be extended to all 1D cases.
Catherine Rio

File size: 32.3 KB
Line 
1!----------------------------------------------------------------------
2! forcing_les = .T. : Impose a constant cooling
3! forcing_radconv = .T. : Pure radiative-convective equilibrium:
4!----------------------------------------------------------------------
5
6
7      nq1=0
8      nq2=0
9
10      if (forcing_les .or. forcing_radconv                                      &
11     &    .or. forcing_GCSSold .or. forcing_fire) then
12
13      if (forcing_fire) then
14!----------------------------------------------------------------------
15!read fire forcings from fire.nc
16!----------------------------------------------------------------------
17      fich_fire='fire.nc'
18      call read_fire(fich_fire,nlev_fire,nt_fire                                &
19     &     ,height,tttprof,qtprof,uprof,vprof,e12prof                           &
20     &     ,ugprof,vgprof,wfls,dqtdxls                                          &
21     &     ,dqtdyls,dqtdtls,thlpcar)
22      write(*,*) 'Forcing FIRE lu'
23      kmax=120            ! nombre de niveaux dans les profils et forcages
24      else 
25!----------------------------------------------------------------------
26! Read profiles from files: prof.inp.001 and lscale.inp.001
27! (repris de readlesfiles)
28!----------------------------------------------------------------------
29
30      call readprofiles(nlev_max,kmax,nqtot,height,                             &
31     &           tttprof,qtprof,uprof,vprof,                                    &
32     &           e12prof,ugprof,vgprof,                                         &
33     &           wfls,dqtdxls,dqtdyls,dqtdtls,                                  &
34     &           thlpcar,qprof,nq1,nq2)
35      endif
36
37! compute altitudes of play levels.
38      zlay(1) =zsurf +  rd*tsurf*(psurf-play(1))/(rg*psurf)
39      do l = 2,llm
40        zlay(l) = zlay(l-1)+rd*tsurf*(psurf-play(1))/(rg*psurf)
41      enddo
42
43!----------------------------------------------------------------------
44! Interpolation of the profiles given on the input file to
45! model levels
46!----------------------------------------------------------------------
47      zlay(1) = zsurf +  rd*tsurf*(psurf-play(1))/(rg*psurf)
48      do l=1,llm
49        ! Above the max altutide of the input file
50
51        if (zlay(l)<height(kmax)) mxcalc=l
52
53        frac = (height(kmax)-zlay(l))/(height (kmax)-height(kmax-1))
54        ttt =tttprof(kmax)-frac*(tttprof(kmax)-tttprof(kmax-1))
55       if ((forcing_GCSSold .AND. tp_ini_GCSSold) .OR. forcing_fire)then ! pot. temp. in initial profile
56          temp(l) = ttt*(play(l)/pzero)**rkappa
57          teta(l) = ttt
58       else
59          temp(l) = ttt
60          teta(l) = ttt*(pzero/play(l))**rkappa
61       endif
62          print *,' temp,teta ',l,temp(l),teta(l)
63        q(l,1)  = qtprof(kmax)-frac*( qtprof(kmax)- qtprof(kmax-1))
64        u(l)    =  uprof(kmax)-frac*(  uprof(kmax)-  uprof(kmax-1))
65        v(l)    =  vprof(kmax)-frac*(  vprof(kmax)-  vprof(kmax-1))
66        ug(l)   = ugprof(kmax)-frac*( ugprof(kmax)- ugprof(kmax-1))
67        vg(l)   = vgprof(kmax)-frac*( vgprof(kmax)- vgprof(kmax-1))
68        IF (nq2>0) q(l,nq1:nq2)=qprof(kmax,nq1:nq2)                         &
69     &               -frac*(qprof(kmax,nq1:nq2)-qprof(kmax-1,nq1:nq2))
70        omega(l)=   wfls(kmax)-frac*(   wfls(kmax)-   wfls(kmax-1))
71
72        dq_dyn(l,1) = dqtdtls(kmax)-frac*(dqtdtls(kmax)-dqtdtls(kmax-1))
73        dt_cooling(l)=thlpcar(kmax)-frac*(thlpcar(kmax)-thlpcar(kmax-1))
74        do k=2,kmax
75          frac = (height(k)-zlay(l))/(height(k)-height(k-1))
76          if(l==1) print*,'k, height, tttprof',k,height(k),tttprof(k)
77          if(zlay(l)>height(k-1).and.zlay(l)<height(k)) then
78            ttt =tttprof(k)-frac*(tttprof(k)-tttprof(k-1))
79       if ((forcing_GCSSold .AND. tp_ini_GCSSold) .OR. forcing_fire)then ! pot. temp. in initial profile
80          temp(l) = ttt*(play(l)/pzero)**rkappa
81          teta(l) = ttt
82       else
83          temp(l) = ttt
84          teta(l) = ttt*(pzero/play(l))**rkappa
85       endif
86          print *,' temp,teta ',l,temp(l),teta(l)
87            q(l,1)  = qtprof(k)-frac*( qtprof(k)- qtprof(k-1))
88            u(l)    =  uprof(k)-frac*(  uprof(k)-  uprof(k-1))
89            v(l)    =  vprof(k)-frac*(  vprof(k)-  vprof(k-1))
90            ug(l)   = ugprof(k)-frac*( ugprof(k)- ugprof(k-1))
91            vg(l)   = vgprof(k)-frac*( vgprof(k)- vgprof(k-1))
92            IF (nq2>0) q(l,nq1:nq2)=qprof(k,nq1:nq2)                        &
93     &                   -frac*(qprof(k,nq1:nq2)-qprof(k-1,nq1:nq2))
94            omega(l)=   wfls(k)-frac*(   wfls(k)-   wfls(k-1))
95            dq_dyn(l,1)=dqtdtls(k)-frac*(dqtdtls(k)-dqtdtls(k-1))
96            dt_cooling(l)=thlpcar(k)-frac*(thlpcar(k)-thlpcar(k-1))
97          elseif(zlay(l)<height(1)) then ! profils uniformes pour z<height(1)
98            ttt =tttprof(1)
99       if ((forcing_GCSSold .AND. tp_ini_GCSSold) .OR. forcing_fire)then ! pot. temp. in initial profile
100          temp(l) = ttt*(play(l)/pzero)**rkappa
101          teta(l) = ttt
102       else
103          temp(l) = ttt
104          teta(l) = ttt*(pzero/play(l))**rkappa
105       endif
106            q(l,1)  = qtprof(1)
107            u(l)    =  uprof(1)
108            v(l)    =  vprof(1)
109            ug(l)   = ugprof(1)
110            vg(l)   = vgprof(1)
111            omega(l)=   wfls(1)
112            IF (nq2>0) q(l,nq1:nq2)=qprof(1,nq1:nq2)
113            dq_dyn(l,1)  =dqtdtls(1)
114            dt_cooling(l)=thlpcar(1)
115          endif
116        enddo
117
118        temp(l)=max(min(temp(l),350.),150.)
119        rho(l)  = play(l)/(rd*temp(l)*(1.+(rv/rd-1.)*q(l,1)))
120        if (l .lt. llm) then
121          zlay(l+1) = zlay(l) + (play(l)-play(l+1))/(rg*rho(l))
122        endif
123        omega2(l)=-rho(l)*omega(l)
124        omega(l)= omega(l)*(-rg*rho(l)) !en Pa/s
125        if (l>1) then
126        if(zlay(l-1)>height(kmax)) then
127           omega(l)=0.0
128           omega2(l)=0.0
129        endif   
130        endif
131        if(q(l,1)<0.) q(l,1)=0.0
132        q(l,2)  = 0.0
133      enddo
134
135      endif ! forcing_les .or. forcing_GCSSold .or. forcing_fire
136!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
137!---------------------------------------------------------------------
138! Forcing for GCSSold:
139!---------------------------------------------------------------------
140      if (forcing_GCSSold) then
141       fich_gcssold_ctl = './forcing.ctl'
142       fich_gcssold_dat = './forcing8.dat'
143       call copie(llm,play,psurf,fich_gcssold_ctl)
144       call get_uvd2(it,timestep,fich_gcssold_ctl,fich_gcssold_dat,         &
145     &               ht_gcssold,hq_gcssold,hw_gcssold,                      &
146     &               hu_gcssold,hv_gcssold,                                 &
147     &               hthturb_gcssold,hqturb_gcssold,Ts_gcssold,             &
148     &               imp_fcg_gcssold,ts_fcg_gcssold,                        &
149     &               Tp_fcg_gcssold,Turb_fcg_gcssold)
150       print *,' get_uvd2 -> hqturb_gcssold ',hqturb_gcssold
151      endif ! forcing_GCSSold
152!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
153!---------------------------------------------------------------------
154! Forcing for RICO:
155!---------------------------------------------------------------------
156      if (forcing_rico) then
157
158!       call writefield_phy('omega', omega,llm+1)
159      fich_rico = 'rico.txt'
160       call read_rico(fich_rico,nlev_rico,ps_rico,play                      &
161     &             ,ts_rico,t_rico,q_rico,u_rico,v_rico,w_rico              &
162     &             ,dth_rico,dqh_rico)
163        print*, ' on a lu et prepare RICO'
164
165       mxcalc=llm
166       print *, airefi, ' airefi '
167       do l = 1, llm
168       rho(l)  = play(l)/(rd*t_rico(l)*(1.+(rv/rd-1.)*q_rico(l)))
169       temp(l) = t_rico(l)
170       q(l,1) = q_rico(l)
171       q(l,2) = 0.0
172       u(l) = u_rico(l)
173       v(l) = v_rico(l)
174       ug(l)=u_rico(l)
175       vg(l)=v_rico(l)
176       omega(l) = -w_rico(l)*rg
177       omega2(l) = omega(l)/rg*airefi
178       enddo
179      endif
180!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
181!---------------------------------------------------------------------
182! Forcing from TOGA-COARE experiment (Ciesielski et al. 2002) :
183!---------------------------------------------------------------------
184
185      if (forcing_toga) then
186
187! read TOGA-COARE forcing (native vertical grid, nt_toga timesteps):
188      fich_toga = './d_toga/ifa_toga_coare_v21_dime.txt'
189      CALL read_togacoare(fich_toga,nlev_toga,nt_toga                       &
190     &         ,ts_toga,plev_toga,t_toga,q_toga,u_toga,v_toga,w_toga        &
191     &         ,ht_toga,vt_toga,hq_toga,vq_toga)
192
193       write(*,*) 'Forcing TOGA lu'
194
195! time interpolation for initial conditions:
196      write(*,*) 'AVT 1ere INTERPOLATION: day,day1 = ',day,day1
197      CALL interp_toga_time(daytime,day1,annee_ref                          &
198     &             ,year_ini_toga,day_ju_ini_toga,nt_toga,dt_toga           &
199     &             ,nlev_toga,ts_toga,plev_toga,t_toga,q_toga,u_toga        &
200     &             ,v_toga,w_toga,ht_toga,vt_toga,hq_toga,vq_toga           &
201     &             ,ts_prof,plev_prof,t_prof,q_prof,u_prof,v_prof,w_prof    &
202     &             ,ht_prof,vt_prof,hq_prof,vq_prof)
203
204! vertical interpolation:
205      CALL interp_toga_vertical(play,nlev_toga,plev_prof                    &
206     &         ,t_prof,q_prof,u_prof,v_prof,w_prof                          &
207     &         ,ht_prof,vt_prof,hq_prof,vq_prof                             &
208     &         ,t_mod,q_mod,u_mod,v_mod,w_mod                               &
209     &         ,ht_mod,vt_mod,hq_mod,vq_mod,mxcalc)
210       write(*,*) 'Profil initial forcing TOGA interpole'
211
212! initial and boundary conditions :
213      tsurf = ts_prof
214      write(*,*) 'SST initiale: ',tsurf
215      do l = 1, llm
216       temp(l) = t_mod(l)
217       q(l,1) = q_mod(l)
218       q(l,2) = 0.0
219       u(l) = u_mod(l)
220       v(l) = v_mod(l)
221       omega(l) = w_mod(l)
222       omega2(l)=omega(l)/rg*airefi ! flxmass_w calcule comme ds physiq
223!?       rho(l)  = play(l)/(rd*temp(l)*(1.+(rv/rd-1.)*q(l,1)))
224!?       omega2(l)=-rho(l)*omega(l)
225       alpha = rd*temp(l)*(1.+(rv/rd-1.)*q(l,1))/play(l)
226       d_th_adv(l) = alpha*omega(l)/rcpd-(ht_mod(l)+vt_mod(l))
227       d_q_adv(l,1) = -(hq_mod(l)+vq_mod(l))
228       d_q_adv(l,2) = 0.0
229      enddo
230
231      endif ! forcing_toga
232!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
233!---------------------------------------------------------------------
234! Forcing from TWPICE experiment (Shaocheng et al. 2010) :
235!---------------------------------------------------------------------
236
237      if (forcing_twpice) then
238!read TWP-ICE forcings
239     fich_twpice='d_twpi/twp180iopsndgvarana_v2.1_C3.c1.20060117.000000.cdf'
240      call read_twpice(fich_twpice,nlev_twpi,nt_twpi                        &
241     &     ,ts_twpi,plev_twpi,t_twpi,q_twpi,u_twpi,v_twpi,w_twpi            &
242     &     ,ht_twpi,vt_twpi,hq_twpi,vq_twpi)
243
244      write(*,*) 'Forcing TWP-ICE lu'
245!Time interpolation for initial conditions using TOGA interpolation routine
246         write(*,*) 'AVT 1ere INTERPOLATION: day,day1 = ',daytime,day1   
247      CALL interp_toga_time(daytime,day1,annee_ref                          &
248     &          ,year_ini_twpi,day_ju_ini_twpi,nt_twpi,dt_twpi,nlev_twpi    &
249     &             ,ts_twpi,plev_twpi,t_twpi,q_twpi,u_twpi,v_twpi,w_twpi    &
250     &             ,ht_twpi,vt_twpi,hq_twpi,vq_twpi                         &
251     &             ,ts_proftwp,plev_proftwp,t_proftwp,q_proftwp             &
252     &             ,u_proftwp,v_proftwp,w_proftwp                           &
253     &             ,ht_proftwp,vt_proftwp,hq_proftwp,vq_proftwp)
254
255! vertical interpolation using TOGA interpolation routine:
256!      write(*,*)'avant interp vert', t_proftwp
257      CALL interp_toga_vertical(play,nlev_twpi,plev_proftwp                 &
258     &         ,t_proftwp,q_proftwp,u_proftwp,v_proftwp,w_proftwp           &
259     &         ,ht_proftwp,vt_proftwp,hq_proftwp,vq_proftwp                 &
260     &         ,t_mod,q_mod,u_mod,v_mod,w_mod                               &
261     &         ,ht_mod,vt_mod,hq_mod,vq_mod,mxcalc)
262!       write(*,*) 'Profil initial forcing TWP-ICE interpole',t_mod
263
264! initial and boundary conditions :
265!      tsurf = ts_proftwp
266      write(*,*) 'SST initiale: ',tsurf
267      do l = 1, llm
268       temp(l) = t_mod(l)
269       q(l,1) = q_mod(l)
270       q(l,2) = 0.0
271       u(l) = u_mod(l)
272       v(l) = v_mod(l)
273       omega(l) = w_mod(l)
274       omega2(l)=omega(l)/rg*airefi ! flxmass_w calcule comme ds physiq
275
276       alpha = rd*temp(l)*(1.+(rv/rd-1.)*q(l,1))/play(l)
277!on applique le forcage total au premier pas de temps
278!attention: signe different de toga
279       d_th_adv(l) = alpha*omega(l)/rcpd+(ht_mod(l)+vt_mod(l))
280       d_q_adv(l,1) = (hq_mod(l)+vq_mod(l))
281       d_q_adv(l,2) = 0.0
282      enddo     
283       
284      endif !forcing_twpice
285
286!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
287!---------------------------------------------------------------------
288! Forcing from AMMA experiment (Couvreux et al. 2010) :
289!---------------------------------------------------------------------
290
291      if (forcing_amma) then
292
293      call read_1D_cases
294
295      write(*,*) 'Forcing AMMA lu'
296
297!champs initiaux:
298      do k=1,nlev_amma
299         th_ammai(k)=th_amma(k)
300         q_ammai(k)=q_amma(k)
301         u_ammai(k)=u_amma(k)
302         v_ammai(k)=v_amma(k)
303         vitw_ammai(k)=vitw_amma(k,12)
304         ht_ammai(k)=ht_amma(k,12)
305         hq_ammai(k)=hq_amma(k,12)
306         vt_ammai(k)=0.
307         vq_ammai(k)=0.
308      enddo   
309      omega(:)=0.     
310      omega2(:)=0.
311      rho(:)=0.
312! vertical interpolation using TOGA interpolation routine:
313!      write(*,*)'avant interp vert', t_proftwp
314      CALL interp_toga_vertical(play,nlev_amma,plev_amma                    &
315     &         ,th_ammai,q_ammai,u_ammai,v_ammai,vitw_ammai                 &
316     &         ,ht_ammai,vt_ammai,hq_ammai,vq_ammai                         &
317     &         ,t_mod,q_mod,u_mod,v_mod,w_mod                               &
318     &         ,ht_mod,vt_mod,hq_mod,vq_mod,mxcalc)
319!       write(*,*) 'Profil initial forcing TWP-ICE interpole',t_mod
320
321! initial and boundary conditions :
322!      tsurf = ts_proftwp
323      write(*,*) 'SST initiale mxcalc: ',tsurf,mxcalc
324      do l = 1, llm
325! Ligne du dessous à decommenter si on lit theta au lieu de temp
326!      temp(l) = t_mod(l)*(play(l)/pzero)**rkappa
327       temp(l) = t_mod(l) 
328       q(l,1) = q_mod(l)
329       q(l,2) = 0.0
330!      print *,'read_forc: l,temp,q=',l,temp(l),q(l,1)
331       u(l) = u_mod(l)
332       v(l) = v_mod(l)
333       rho(l)  = play(l)/(rd*temp(l)*(1.+(rv/rd-1.)*q(l,1)))
334       omega(l) = w_mod(l)*(-rg*rho(l))
335       omega2(l)=omega(l)/rg*airefi ! flxmass_w calcule comme ds physiq
336
337       alpha = rd*temp(l)*(1.+(rv/rd-1.)*q(l,1))/play(l)
338!on applique le forcage total au premier pas de temps
339!attention: signe different de toga
340       d_th_adv(l) = alpha*omega(l)/rcpd+ht_mod(l)
341!forcage en th
342!       d_th_adv(l) = ht_mod(l)
343       d_q_adv(l,1) = hq_mod(l)
344       d_q_adv(l,2) = 0.0
345       dt_cooling(l)=0.
346      enddo     
347       write(*,*) 'Prof initeforcing AMMA interpole temp39',temp(39)
348     
349
350!     ok_flux_surf=.false.
351      fsens=-1.*sens_amma(12)
352      flat=-1.*lat_amma(12)
353       
354      endif !forcing_amma
355
356
357!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
358!---------------------------------------------------------------------
359! Forcing from DICE experiment (see file DICE_protocol_vn2-3.pdf)
360!---------------------------------------------------------------------
361
362      if (forcing_dice) then
363!read DICE forcings
364      fich_dice='dice_driver.nc'
365      call read_dice(fich_dice,nlev_dice,nt_dice                    &
366     &     ,zz_dice,plev_dice,th_dice,qv_dice,u_dice,v_dice,o3_dice &
367     &     ,shf_dice,lhf_dice,lwup_dice,swup_dice,tg_dice,ustar_dice& 
368     &     ,psurf_dice,ug_dice,vg_dice,ht_dice,hq_dice              &
369     &     ,hu_dice,hv_dice,w_dice,omega_dice)
370
371      write(*,*) 'Forcing DICE lu'
372
373!champs initiaux:
374      do k=1,nlev_dice
375         th_dicei(k)=th_dice(k)
376         qv_dicei(k)=qv_dice(k)
377         u_dicei(k)=u_dice(k)
378         v_dicei(k)=v_dice(k)
379         o3_dicei(k)=o3_dice(k)
380         ht_dicei(k)=ht_dice(k,1)
381         hq_dicei(k)=hq_dice(k,1)
382         hu_dicei(k)=hu_dice(k,1)
383         hv_dicei(k)=hv_dice(k,1)
384         w_dicei(k)=w_dice(k,1)
385         omega_dicei(k)=omega_dice(k,1)
386      enddo   
387      omega(:)=0.     
388      omega2(:)=0.
389      rho(:)=0.
390! vertical interpolation using TOGA interpolation routine:
391!      write(*,*)'avant interp vert', t_proftwp
392!
393!     CALL interp_dice_time(daytime,day1,annee_ref
394!    i             ,year_ini_dice,day_ju_ini_dice,nt_dice,dt_dice
395!    i             ,nlev_dice,shf_dice,lhf_dice,lwup_dice,swup_dice
396!    i             ,tg_dice,ustar_dice,psurf_dice,ug_dice,vg_dice
397!    i             ,ht_dice,hq_dice,hu_dice,hv_dice,w_dice,omega_dice
398!    o             ,shf_prof,lhf_prof,lwup_prof,swup_prof,tg_prof
399!    o             ,ustar_prof,psurf_prof,ug_profd,vg_profd
400!    o             ,ht_profd,hq_profd,hu_profd,hv_profd,w_profd
401!    o             ,omega_profd)
402
403      CALL interp_dice_vertical(play,nlev_dice,nt_dice,plev_dice       &
404     &         ,th_dicei,qv_dicei,u_dicei,v_dicei,o3_dicei             &
405     &         ,ht_dicei,hq_dicei,hu_dicei,hv_dicei,w_dicei,omega_dicei&
406     &         ,th_mod,qv_mod,u_mod,v_mod,o3_mod                       &
407     &         ,ht_mod,hq_mod,hu_mod,hv_mod,w_mod,omega_mod,mxcalc)
408
409! Pour tester les advections horizontales de T et Q, on met w_mod et omega_mod à zero (MPL 20131108)
410!     w_mod(:,:)=0.
411!     omega_mod(:,:)=0.
412
413!       write(*,*) 'Profil initial forcing DICE interpole',t_mod
414! Les forcages DICE sont donnes /jour et non /seconde !
415      ht_mod(:)=ht_mod(:)/86400.
416      hq_mod(:)=hq_mod(:)/86400.
417      hu_mod(:)=hu_mod(:)/86400.
418      hv_mod(:)=hv_mod(:)/86400.
419
420! initial and boundary conditions :
421      write(*,*) 'SST initiale mxcalc: ',tsurf,mxcalc
422      do l = 1, llm
423! Ligne du dessous à decommenter si on lit theta au lieu de temp
424       temp(l) = th_mod(l)*(play(l)/pzero)**rkappa
425!      temp(l) = t_mod(l) 
426       q(l,1) = qv_mod(l)
427       q(l,2) = 0.0
428!      print *,'read_forc: l,temp,q=',l,temp(l),q(l,1)
429       u(l) = u_mod(l)
430       v(l) = v_mod(l)
431       ug(l)=ug_dice(1)
432       vg(l)=vg_dice(1)
433       rho(l)  = play(l)/(rd*temp(l)*(1.+(rv/rd-1.)*q(l,1)))
434!      omega(l) = w_mod(l)*(-rg*rho(l))
435       omega(l) = omega_mod(l)
436       omega2(l)=omega(l)/rg*airefi ! flxmass_w calcule comme ds physiq
437
438       alpha = rd*temp(l)*(1.+(rv/rd-1.)*q(l,1))/play(l)
439!on applique le forcage total au premier pas de temps
440!attention: signe different de toga
441       d_th_adv(l) = alpha*omega(l)/rcpd+ht_mod(l)
442!forcage en th
443!       d_th_adv(l) = ht_mod(l)
444       d_q_adv(l,1) = hq_mod(l)
445       d_q_adv(l,2) = 0.0
446       dt_cooling(l)=0.
447      enddo     
448       write(*,*) 'Profil initial forcing DICE interpole temp39',temp(39)
449     
450
451!     ok_flux_surf=.false.
452      fsens=-1.*shf_dice(1)
453      flat=-1.*lhf_dice(1)
454! Le cas Dice doit etre force avec ustar mais on peut simplifier en forcant par
455! le coefficient de trainee en surface cd**2=ustar*vent(k=1)
456! On commence ici a stocker ustar dans cdrag puis on terminera le calcul dans pbl_surface
457! MPL 05082013
458      ust=ustar_dice(1)
459      tg=tg_dice(1)
460      print *,'ust= ',ust
461      IF (tsurf .LE. 0.) THEN
462       tsurf= tg_dice(1)
463      ENDIF
464      psurf= psurf_dice(1)
465      solsw_in = (1.-albedo)/albedo*swup_dice(1)
466      sollw_in = (0.7*RSIGMA*temp(1)**4)-lwup_dice(1)
467      PRINT *,'1D_READ_FORC : solsw, sollw',solsw_in,sollw_in
468      endif !forcing_dice
469
470!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
471!---------------------------------------------------------------------
472! Forcing from Arm_Cu case                   
473! For this case, ifa_armcu.txt contains sensible, latent heat fluxes
474! large scale advective forcing,radiative forcing
475! and advective tendency of theta and qt to be applied
476!---------------------------------------------------------------------
477
478      if (forcing_armcu) then
479! read armcu forcing :
480       write(*,*) 'Avant lecture Forcing Arm_Cu'
481      fich_armcu = './ifa_armcu.txt'
482      CALL read_armcu(fich_armcu,nlev_armcu,nt_armcu,                       &
483     & sens_armcu,flat_armcu,adv_theta_armcu,                               &
484     & rad_theta_armcu,adv_qt_armcu)
485       write(*,*) 'Forcing Arm_Cu lu'
486
487!----------------------------------------------------------------------
488! Read profiles from file: prof.inp.19 or prof.inp.40 
489! For this case, profiles are given for two vertical resolution
490! 19 or 40 levels
491!
492! Comment from: http://www.knmi.nl/samenw/eurocs/ARM/profiles.html
493! Note that the initial profiles contain no liquid water! 
494! (so potential temperature can be interpreted as liquid water
495! potential temperature and water vapor as total water)
496! profiles are given at full levels
497!----------------------------------------------------------------------
498
499      call readprofile_armcu(nlev_max,kmax,height,play_mod,u_mod,           &
500     &           v_mod,theta_mod,t_mod,qv_mod,rv_mod,ap,bp)
501
502! time interpolation for initial conditions:
503      write(*,*) 'AVT 1ere INTERPOLATION: day,day1 = ',day,day1
504
505      print *,'Avant interp_armcu_time'
506      print *,'daytime=',daytime
507      print *,'day1=',day1
508      print *,'annee_ref=',annee_ref
509      print *,'year_ini_armcu=',year_ini_armcu
510      print *,'day_ju_ini_armcu=',day_ju_ini_armcu
511      print *,'nt_armcu=',nt_armcu
512      print *,'dt_armcu=',dt_armcu
513      print *,'nlev_armcu=',nlev_armcu
514      CALL interp_armcu_time(daytime,day1,annee_ref                         &
515     &            ,year_ini_armcu,day_ju_ini_armcu,nt_armcu,dt_armcu        &
516     &            ,nlev_armcu,sens_armcu,flat_armcu,adv_theta_armcu         &
517     &            ,rad_theta_armcu,adv_qt_armcu,sens_prof,flat_prof         &
518     &            ,adv_theta_prof,rad_theta_prof,adv_qt_prof)
519       write(*,*) 'Forcages interpoles dans temps'
520
521! No vertical interpolation if nlev imposed to 19 or 40
522! The vertical grid stops at 4000m # 600hPa
523      mxcalc=llm
524
525! initial and boundary conditions :
526!     tsurf = ts_prof
527! tsurf read in lmdz1d.def
528      write(*,*) 'Tsurf initiale: ',tsurf
529      do l = 1, llm
530       play(l)=play_mod(l)*100.
531       presnivs(l)=play(l)
532       zlay(l)=height(l)
533       temp(l) = t_mod(l)
534       teta(l)=theta_mod(l)
535       q(l,1) = qv_mod(l)/1000.
536! No liquid water in the initial profil
537       q(l,2) = 0.
538       u(l) = u_mod(l)
539       ug(l)= u_mod(l)
540       v(l) = v_mod(l)
541       vg(l)= v_mod(l)
542! Advective forcings are given in K or g/kg ... per HOUR
543!      IF(height(l).LT.1000) THEN
544!        d_th_adv(l) = (adv_theta_prof + rad_theta_prof)/3600.
545!        d_q_adv(l,1) = adv_qt_prof/1000./3600.
546!        d_q_adv(l,2) = 0.0
547!      ELSEIF (height(l).GE.1000.AND.height(l).LT.3000) THEN
548!        d_th_adv(l) = (adv_theta_prof + rad_theta_prof)*
549!    :               (1-(height(l)-1000.)/2000.)
550!        d_th_adv(l) = d_th_adv(l)/3600.
551!        d_q_adv(l,1) = adv_qt_prof*(1-(height(l)-1000.)/2000.)
552!        d_q_adv(l,1) = d_q_adv(l,1)/1000./3600.
553!        d_q_adv(l,2) = 0.0
554!      ELSE
555!        d_th_adv(l) = 0.0
556!        d_q_adv(l,1) = 0.0
557!        d_q_adv(l,2) = 0.0
558!      ENDIF
559      enddo
560! plev at half levels is given in proh.inp.19 or proh.inp.40 files
561      plev(1)= ap(llm+1)+bp(llm+1)*psurf
562      do l = 1, llm
563      plev(l+1) = ap(llm-l+1)+bp(llm-l+1)*psurf
564      print *,'Read_forc: l height play plev zlay temp',                    &
565     &   l,height(l),play(l),plev(l),zlay(l),temp(l)
566      enddo
567! For this case, fluxes are imposed
568       fsens=-1*sens_prof
569       flat=-1*flat_prof
570
571      endif ! forcing_armcu
572!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
573!---------------------------------------------------------------------
574! Forcing from transition case of Irina Sandu                 
575!---------------------------------------------------------------------
576
577      if (forcing_sandu) then
578       write(*,*) 'Avant lecture Forcing SANDU'
579
580! read sanduref forcing :
581      fich_sandu = './ifa_sanduref.txt'
582      CALL read_sandu(fich_sandu,nlev_sandu,nt_sandu,ts_sandu)
583
584       write(*,*) 'Forcing SANDU lu'
585
586!----------------------------------------------------------------------
587! Read profiles from file: prof.inp.001 
588!----------------------------------------------------------------------
589
590      call readprofile_sandu(nlev_max,kmax,height,plev_profs,t_profs,       &
591     &           thl_profs,q_profs,u_profs,v_profs,                         &
592     &           w_profs,omega_profs,o3mmr_profs)
593
594! time interpolation for initial conditions:
595      write(*,*) 'AVT 1ere INTERPOLATION: day,day1 = ',day,day1
596! ATTENTION, cet appel ne convient pas pour le cas SANDU !!
597! revoir 1DUTILS.h et les arguments
598
599      print *,'Avant interp_sandu_time'
600      print *,'daytime=',daytime
601      print *,'day1=',day1
602      print *,'annee_ref=',annee_ref
603      print *,'year_ini_sandu=',year_ini_sandu
604      print *,'day_ju_ini_sandu=',day_ju_ini_sandu
605      print *,'nt_sandu=',nt_sandu
606      print *,'dt_sandu=',dt_sandu
607      print *,'nlev_sandu=',nlev_sandu
608      CALL interp_sandu_time(daytime,day1,annee_ref                         &
609     &             ,year_ini_sandu,day_ju_ini_sandu,nt_sandu,dt_sandu       &
610     &             ,nlev_sandu                                              &
611     &             ,ts_sandu,ts_prof)
612
613! vertical interpolation:
614      print *,'Avant interp_vertical: nlev_sandu=',nlev_sandu
615      CALL interp_sandu_vertical(play,nlev_sandu,plev_profs                 &
616     &         ,t_profs,thl_profs,q_profs,u_profs,v_profs,w_profs           &
617     &         ,omega_profs,o3mmr_profs                                     &
618     &         ,t_mod,thl_mod,q_mod,u_mod,v_mod,w_mod                       &
619     &         ,omega_mod,o3mmr_mod,mxcalc)
620       write(*,*) 'Profil initial forcing SANDU interpole'
621
622! initial and boundary conditions :
623      tsurf = ts_prof
624      write(*,*) 'SST initiale: ',tsurf
625      do l = 1, llm
626       temp(l) = t_mod(l)
627       tetal(l)=thl_mod(l)
628       q(l,1) = q_mod(l)
629       q(l,2) = 0.0
630       u(l) = u_mod(l)
631       v(l) = v_mod(l)
632       w(l) = w_mod(l)
633       omega(l) = omega_mod(l)
634       omega2(l)=omega(l)/rg*airefi ! flxmass_w calcule comme ds physiq
635!?       rho(l)  = play(l)/(rd*temp(l)*(1.+(rv/rd-1.)*q(l,1)))
636!?       omega2(l)=-rho(l)*omega(l)
637       alpha = rd*temp(l)*(1.+(rv/rd-1.)*q(l,1))/play(l)
638!      d_th_adv(l) = alpha*omega(l)/rcpd+vt_mod(l)
639!      d_q_adv(l,1) = vq_mod(l)
640       d_th_adv(l) = alpha*omega(l)/rcpd
641       d_q_adv(l,1) = 0.0
642       d_q_adv(l,2) = 0.0
643      enddo
644
645      endif ! forcing_sandu
646!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
647!---------------------------------------------------------------------
648! Forcing from Astex case
649!---------------------------------------------------------------------
650
651      if (forcing_astex) then
652       write(*,*) 'Avant lecture Forcing Astex'
653
654! read astex forcing :
655      fich_astex = './ifa_astex.txt'
656      CALL read_astex(fich_astex,nlev_astex,nt_astex,div_astex,ts_astex,    &
657     &  ug_astex,vg_astex,ufa_astex,vfa_astex)
658
659       write(*,*) 'Forcing Astex lu'
660
661!----------------------------------------------------------------------
662! Read profiles from file: prof.inp.001
663!----------------------------------------------------------------------
664
665      call readprofile_astex(nlev_max,kmax,height,plev_profa,t_profa,       &
666     &           thl_profa,qv_profa,ql_profa,qt_profa,u_profa,v_profa,      &
667     &           w_profa,tke_profa,o3mmr_profa)
668
669! time interpolation for initial conditions:
670      write(*,*) 'AVT 1ere INTERPOLATION: day,day1 = ',day,day1
671! ATTENTION, cet appel ne convient pas pour le cas SANDU !!
672! revoir 1DUTILS.h et les arguments
673
674      print *,'Avant interp_astex_time'
675      print *,'daytime=',daytime
676      print *,'day1=',day1
677      print *,'annee_ref=',annee_ref
678      print *,'year_ini_astex=',year_ini_astex
679      print *,'day_ju_ini_astex=',day_ju_ini_astex
680      print *,'nt_astex=',nt_astex
681      print *,'dt_astex=',dt_astex
682      print *,'nlev_astex=',nlev_astex
683      CALL interp_astex_time(daytime,day1,annee_ref                         &
684     &             ,year_ini_astex,day_ju_ini_astex,nt_astex,dt_astex       &
685     &             ,nlev_astex,div_astex,ts_astex,ug_astex,vg_astex         &
686     &             ,ufa_astex,vfa_astex,div_prof,ts_prof,ug_prof,vg_prof    &
687     &             ,ufa_prof,vfa_prof)
688
689! vertical interpolation:
690      print *,'Avant interp_vertical: nlev_astex=',nlev_astex
691      CALL interp_astex_vertical(play,nlev_astex,plev_profa                 &
692     &         ,t_profa,thl_profa,qv_profa,ql_profa,qt_profa                &
693     &         ,u_profa,v_profa,w_profa,tke_profa,o3mmr_profa               &
694     &         ,t_mod,thl_mod,qv_mod,ql_mod,qt_mod,u_mod,v_mod,w_mod        &
695     &         ,tke_mod,o3mmr_mod,mxcalc)
696       write(*,*) 'Profil initial forcing Astex interpole'
697
698! initial and boundary conditions :
699      tsurf = ts_prof
700      write(*,*) 'SST initiale: ',tsurf
701      do l = 1, llm
702       temp(l) = t_mod(l)
703       tetal(l)=thl_mod(l)
704       q(l,1) = qv_mod(l)
705       q(l,2) = ql_mod(l)
706       u(l) = u_mod(l)
707       v(l) = v_mod(l)
708       w(l) = w_mod(l)
709       omega(l) = w_mod(l)
710!      omega2(l)=omega(l)/rg*airefi ! flxmass_w calcule comme ds physiq
711!      rho(l)  = play(l)/(rd*temp(l)*(1.+(rv/rd-1.)*q(l,1)))
712!      omega2(l)=-rho(l)*omega(l)
713       alpha = rd*temp(l)*(1.+(rv/rd-1.)*q(l,1))/play(l)
714!      d_th_adv(l) = alpha*omega(l)/rcpd+vt_mod(l)
715!      d_q_adv(l,1) = vq_mod(l)
716       d_th_adv(l) = alpha*omega(l)/rcpd
717       d_q_adv(l,1) = 0.0
718       d_q_adv(l,2) = 0.0
719      enddo
720
721      endif ! forcing_astex
722!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
723!---------------------------------------------------------------------
724! Forcing from standard case :
725!---------------------------------------------------------------------
726
727      if (forcing_case) then
728
729         write(*,*),'avant call read_1D_cas'
730         call read_1D_cas
731         write(*,*) 'Forcing read'
732
733!Time interpolation for initial conditions using TOGA interpolation routine
734         write(*,*) 'AVT 1ere INTERPOLATION: day,day1 = ',daytime,day1   
735      CALL interp_case_time(day,day1,annee_ref                &
736     &         ,year_ini_cas,day_ju_ini_cas,nt_cas,pdt_cas,nlev_cas       &
737     &         ,ts_cas,plev_cas,t_cas,q_cas,u_cas,v_cas               &
738     &         ,ug_cas,vg_cas,vitw_cas,du_cas,hu_cas,vu_cas           &
739     &         ,dv_cas,hv_cas,vv_cas,dt_cas,ht_cas,vt_cas,dtrad_cas             &
740     &         ,dq_cas,hq_cas,vq_cas,lat_cas,sens_cas                 &
741     &         ,ts_prof_cas,plev_prof_cas,t_prof_cas,q_prof_cas,u_prof_cas,v_prof_cas         &
742     &         ,ug_prof_cas,vg_prof_cas,vitw_prof_cas,du_prof_cas,hu_prof_cas,vu_prof_cas     &
743     &         ,dv_prof_cas,hv_prof_cas,vv_prof_cas,dt_prof_cas,ht_prof_cas,vt_prof_cas,dtrad_prof_cas       &
744     &         ,dq_prof_cas,hq_prof_cas,vq_prof_cas,lat_prof_cas,sens_prof_cas)
745
746! vertical interpolation using TOGA interpolation routine:
747!      write(*,*)'avant interp vert', t_prof
748      CALL interp_case_vertical(play,nlev_cas,plev_prof_cas            &
749     &         ,t_prof_cas,q_prof_cas,u_prof_cas,v_prof_cas,ug_prof_cas,vg_prof_cas,vitw_prof_cas    &
750     &         ,du_prof_cas,hu_prof_cas,vu_prof_cas,dv_prof_cas,hv_prof_cas,vv_prof_cas           &
751     &         ,dt_prof_cas,ht_prof_cas,vt_prof_cas,dtrad_prof_cas,dq_prof_cas,hq_prof_cas,vq_prof_cas           &
752     &         ,t_mod_cas,q_mod_cas,u_mod_cas,v_mod_cas,ug_mod_cas,vg_mod_cas,w_mod_cas           &
753     &         ,du_mod_cas,hu_mod_cas,vu_mod_cas,dv_mod_cas,hv_mod_cas,vv_mod_cas               &
754     &         ,dt_mod_cas,ht_mod_cas,vt_mod_cas,dtrad_mod_cas,dq_mod_cas,hq_mod_cas,vq_mod_cas,mxcalc)
755!       write(*,*) 'Profil initial forcing case interpole',t_mod
756
757! initial and boundary conditions :
758!      tsurf = ts_prof_cas
759      ts_cur = ts_prof_cas
760      psurf=plev_prof_cas(1)
761      write(*,*) 'SST initiale: ',tsurf
762      do l = 1, llm
763       temp(l) = t_mod_cas(l)
764       q(l,1) = q_mod_cas(l)
765       q(l,2) = 0.0
766       u(l) = u_mod_cas(l)
767       v(l) = v_mod_cas(l)
768       omega(l) = w_mod_cas(l)
769       omega2(l)=omega(l)/rg*airefi ! flxmass_w calcule comme ds physiq
770
771       alpha = rd*temp(l)*(1.+(rv/rd-1.)*q(l,1))/play(l)
772!on applique le forcage total au premier pas de temps
773!attention: signe different de toga
774       d_th_adv(l) = alpha*omega(l)/rcpd+(ht_mod_cas(l)+vt_mod_cas(l))
775       d_q_adv(l,1) = (hq_mod_cas(l)+vq_mod_cas(l))
776       d_q_adv(l,2) = 0.0
777       d_u_adv(l) = (hu_mod_cas(l)+vu_mod_cas(l))
778       d_u_adv(l) = (hv_mod_cas(l)+vv_mod_cas(l))
779      enddo     
780       
781      endif !forcing_case
782!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
Note: See TracBrowser for help on using the repository browser.