Ignore:
Timestamp:
Jul 5, 2013, 10:38:13 AM (11 years ago)
Author:
Laurent Fairhead
Message:

Inclusion des cas AMMA et Sandu dans la configuration 1D
MPL


AMMA and Sandu cases included in 1D configuration
MPL

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/phy1d/lmdz1d.F

    r1763 r1780  
    66      use dimphy
    77      use surface_data, only : type_ocean,ok_veget
    8       use pbl_surface_mod, only : pbl_surface_init, pbl_surface_final
     8      use pbl_surface_mod, only : ftsoil, pbl_surface_init,
     9     $                            pbl_surface_final
    910      use fonte_neige_mod, only : fonte_neige_init, fonte_neige_final
    1011
     
    2526#include "compar1d.h"
    2627#include "flux_arp.h"
     28#include "tsoilnudge.h"
    2729#include "fcg_gcssold.h"
    2830!!!#include "fbforcing.h"
     
    8688
    8789        integer :: kmax = llm
    88         integer nlev_max
    89         parameter (nlev_max = 100)
     90        integer nlev_max,llm700
     91        parameter (nlev_max = 1000)
    9092        real timestep, frac, timeit
    9193        real height(nlev_max),tttprof(nlev_max),qtprof(nlev_max),
     
    98100c        integer :: forcing_type
    99101        logical :: forcing_les     = .false.
    100         logical :: forcing_armcu  = .false.
     102        logical :: forcing_armcu   = .false.
    101103        logical :: forcing_rico    = .false.
    102104        logical :: forcing_radconv = .false.
    103105        logical :: forcing_toga    = .false.
    104106        logical :: forcing_twpice  = .false.
     107        logical :: forcing_amma    = .false.
    105108        logical :: forcing_GCM2SCM = .false.
    106109        logical :: forcing_GCSSold = .false.
     110        logical :: forcing_sandu   = .false.
     111        logical :: forcing_astex   = .false.
     112        logical :: forcing_fire    = .false.
    107113        integer :: type_ts_forcing ! 0 = SST constant; 1 = SST read from a file
    108114!                                                            (cf read_tsurf1d.F)
    109115
    110116!vertical advection computation
    111         real d_t_z(llm), d_q_z(llm)
    112         real d_t_dyn_z(llm), d_q_dyn_z(llm)
    113         real zz(llm)
    114         real zfact
     117!       real d_t_z(llm), d_q_z(llm)
     118!       real d_t_dyn_z(llm), d_q_dyn_z(llm)
     119!       real zz(llm)
     120!       real zfact
    115121
    116122!flag forcings
     
    129135      real :: pzero=1.e5
    130136      real :: play (llm),zlay (llm),sig_s(llm),plev(llm+1)
    131       real :: playd(llm),zlayd(llm)
     137      real :: playd(llm),zlayd(llm),ap_amma(llm+1),bp_amma(llm+1),poub
    132138
    133139!---------------------------------------------------------------------
     
    137143      integer :: iq
    138144      real :: phi(llm)
     145      real :: teta(llm),tetal(llm),temp(llm),u(llm),v(llm),w(llm)
    139146      real :: rlat_rad(1),rlon_rad(1)
    140       real :: teta(llm),temp(llm),u(llm),v(llm)
    141147      real :: omega(llm+1),omega2(llm),rho(llm+1)
    142148      real :: ug(llm),vg(llm),fcoriolis
     
    196202!  Fichiers et d'autres variables
    197203!---------------------------------------------------------------------
    198       real ttt
     204      real ttt,bow,q1
    199205      integer :: ierr,k,l,i,it=1,mxcalc
    200206      integer jjmp1
     
    252258!             initial profiles from RICO files
    253259!             LS convergence imposed from RICO files
     260!forcing_type = 6 ==> forcing_amma = .true.
     261!             initial profiles from AMMA nc file
     262!             LS convergence, omega and surface fluxes imposed from AMMA file 
    254263!forcing_type = 40 ==> forcing_GCSSold = .true.
    255264!             initial profile from GCSS file
    256265!             LS convergence imposed from GCSS file
     266!forcing_type = 59 ==> forcing_sandu = .true.
     267!             initial profiles from sanduref file: see prof.inp.001
     268!             SST varying with time and divergence constante: see ifa_sanduref.txt file
     269!             Radiation has to be computed interactively
     270!forcing_type = 60 ==> forcing_astex = .true.
     271!             initial profiles from file: see prof.inp.001
     272!             SST,divergence,ug,vg,ufa,vfa varying with time : see ifa_astex.txt file
     273!             Radiation has to be computed interactively
    257274!forcing_type = 61 ==> forcing_armcu = .true.
    258 !             initial profile from arm_cu file
    259 !             LS convergence imposed from arm_cu file
     275!             initial profiles from file: see prof.inp.001
     276!             sensible and latent heat flux imposed: see ifa_arm_cu_1.txt
     277!             large scale advective forcing & radiative tendencies applied below 1000m: see ifa_arm_cu_2.txt
     278!             use geostrophic wind ug=10m/s vg=0m/s. Duration of the case 53100s
     279!             Radiation to be switched off
    260280!
    261281      if (forcing_type .eq.0) THEN
     
    271291      elseif (forcing_type .eq.5) THEN
    272292       forcing_rico = .true.
     293      elseif (forcing_type .eq.6) THEN
     294       forcing_amma = .true.
    273295      elseif (forcing_type .eq.40) THEN
    274296       forcing_GCSSold = .true.
     297      elseif (forcing_type .eq.59) THEN
     298       forcing_sandu   = .true.
     299      elseif (forcing_type .eq.60) THEN
     300       forcing_astex   = .true.
    275301      elseif (forcing_type .eq.61) THEN
    276302       forcing_armcu = .true.
     
    278304      else
    279305       write (*,*) 'ERROR : unknown forcing_type ', forcing_type
    280        stop 'Forcing_type should be 0,1,2,3 or 40'
     306       stop 'Forcing_type should be 0,1,2,3,4,5,6 or 40,59,60,61'
    281307      ENDIF
    282308      print*,"forcing type=",forcing_type
     
    288314
    289315        type_ts_forcing = 0
    290         if (forcing_toga) type_ts_forcing = 1
     316        if (forcing_toga .or. forcing_sandu .or. forcing_astex)
     317     :    type_ts_forcing = 1
    291318
    292319!---------------------------------------------------------------------
     
    327354c Special case for arm_cu which lasts less than one day : 53100s !! (MPL 20111026)
    328355      IF(forcing_type .EQ. 61) fnday=53100./86400.
     356c Special case for amma which lasts less than one day : 64800s !! (MPL 20120216)
     357      IF(forcing_type .EQ. 6) fnday=64800./86400.
    329358      annee_ref = anneeref
    330359      mois = 1
     
    336365      day_ini = day
    337366      day_end = day_ini + nday
     367
     368      IF (forcing_type .eq.2) THEN
    338369! Convert the initial date of Toga-Coare to Julian day
    339370      call ymds2ju
    340371     $ (year_ini_toga,mth_ini_toga,day_ini_toga,heure,day_ju_ini_toga)
    341372
     373      ELSEIF (forcing_type .eq.4) THEN
    342374! Convert the initial date of TWPICE to Julian day
    343375      call ymds2ju
    344376     $ (year_ini_twpi,mth_ini_twpi,day_ini_twpi,heure_ini_twpi
    345377     $ ,day_ju_ini_twpi)
    346 
    347 ! Convert the initial date of Arm_cu to Julian day
     378      ELSEIF (forcing_type .eq.6) THEN
     379! Convert the initial date of AMMA to Julian day
     380      call ymds2ju
     381     $ (year_ini_amma,mth_ini_amma,day_ini_amma,heure_ini_amma
     382     $ ,day_ju_ini_amma)
     383
     384      ELSEIF (forcing_type .eq.59) THEN
     385! Convert the initial date of Sandu case to Julian day
     386      call ymds2ju
     387     $   (year_ini_sandu,mth_ini_sandu,day_ini_sandu,
     388     $    time_ini*3600.,day_ju_ini_sandu)
     389
     390      ELSEIF (forcing_type .eq.60) THEN
     391! Convert the initial date of Astex case to Julian day
     392      call ymds2ju
     393     $   (year_ini_astex,mth_ini_astex,day_ini_astex,
     394     $    time_ini*3600.,day_ju_ini_astex)
     395
     396      ELSEIF (forcing_type .eq.61) THEN
     397
     398! Convert the initial date of Arm_cu case to Julian day
    348399      call ymds2ju
    349400     $ (year_ini_armcu,mth_ini_armcu,day_ini_armcu,heure_ini_armcu
    350401     $ ,day_ju_ini_armcu)
     402      ENDIF
    351403
    352404      daytime = day + time_ini/24. ! 1st day and initial time of the simulation
     
    435487ccc      zlay=-rd*300.*log(play/psurf)/rg ! moved after reading profiles
    436488
     489      IF (forcing_type .eq. 59) THEN
     490! pour forcing_sandu, on cherche l'indice le plus proche de 700hpa#3000m
    437491      write(*,*) '***********************'
    438492      do l = 1, llm
    439493       write(*,*) 'l,play(l),presnivs(l): ',l,play(l),presnivs(l)
     494       if (trouve_700 .and. play(l).le.70000) then
     495         llm700=l
     496         print *,'llm700,play=',llm700,play(l)/100.
     497         trouve_700= .false.
     498       endif
    440499      enddo
    441500      write(*,*) '***********************'
     501      ENDIF
    442502
    443503c
     
    515575        agesno  = xagesno
    516576        tsoil(:,:,:)=tsurf
     577!------ AMMA 2e run avec modele sol et rayonnement actif (MPL 23052012)
     578!       tsoil(1,1,1)=299.18
     579!       tsoil(1,2,1)=300.08
     580!       tsoil(1,3,1)=301.88
     581!       tsoil(1,4,1)=305.48
     582!       tsoil(1,5,1)=308.00
     583!       tsoil(1,6,1)=308.00
     584!       tsoil(1,7,1)=308.00
     585!       tsoil(1,8,1)=308.00
     586!       tsoil(1,9,1)=308.00
     587!       tsoil(1,10,1)=308.00
     588!       tsoil(1,11,1)=308.00
     589!-----------------------------------------------------------------------
    517590        call pbl_surface_init(qsol, fder, snsrf, qsurfsrf,
    518591     &                                    evap, frugs, agesno, tsoil)
     
    748821       endif
    749822
    750        if (forcing_toga .or. forcing_GCSSold .or. forcing_twpice) then
     823       if (forcing_toga .or. forcing_GCSSold .or. forcing_twpice
     824     :    .or.forcing_amma) then
    751825         fcoriolis=0.0 ; ug=0. ; vg=0.
    752826       endif
     
    813887
    814888        teta=temp*(pzero/play)**rkappa
     889!
     890!---------------------------------------------------------------------
     891!   Nudge soil temperature if requested
     892!---------------------------------------------------------------------
     893
     894      IF (nudge_tsoil) THEN
     895       ftsoil(1,isoil_nudge,:) = ftsoil(1,isoil_nudge,:)
     896     .  -timestep/tau_soil_nudge*(ftsoil(1,isoil_nudge,:)-Tsoil_nudge)
     897      ENDIF
    815898
    816899!---------------------------------------------------------------------
Note: See TracChangeset for help on using the changeset viewer.