MODULE lmdz_scm PRIVATE ! -- We'd love to put IMPLICIT NONE; here... PUBLIC scm CONTAINS SUBROUTINE scm USE ioipsl, ONLY: ju2ymds, ymds2ju, ioconf_calendar, getin USE phys_state_var_mod, ONLY: phys_state_var_init, phys_state_var_end, & clwcon, detr_therm, & qsol, fevap, z0m, z0h, agesno, & du_gwd_rando, du_gwd_front, entr_therm, f0, fm_therm, & falb_dir, falb_dif, & ftsol, beta_aridity, pbl_tke, pctsrf, radsol, rain_fall, snow_fall, ratqs, & rnebcon, rugoro, sig1, w01, solaire_etat0, sollw, sollwdown, & solsw, solswfdiff, t_ancien, q_ancien, u_ancien, v_ancien, rneb_ancien, & wake_delta_pbl_TKE, delta_tsurf, wake_fip, wake_pe, & wake_deltaq, wake_deltat, wake_s, awake_s, wake_dens, & awake_dens, cv_gen, wake_cstar, & zgam, zmax0, zmea, zpic, zsig, & zstd, zthe, zval, ale_bl, ale_bl_trig, alp_bl, ql_ancien, qs_ancien, & prlw_ancien, prsw_ancien, prw_ancien, & u10m, v10m, ale_wake, ale_bl_stat, ratqs_inter_ USE dimphy USE surface_data, ONLY: type_ocean, ok_veget USE pbl_surface_mod, ONLY: ftsoil, pbl_surface_init, & pbl_surface_final USE fonte_neige_mod, ONLY: fonte_neige_init, fonte_neige_final USE infotrac USE control_mod USE indice_sol_mod USE phyaqua_mod USE mod_1D_cases_read_std USE lmdz_print_control, ONLY: lunout, prt_level USE iniphysiq_mod, ONLY: iniphysiq USE mod_const_mpi, ONLY: comm_lmdz USE physiq_mod, ONLY: physiq USE comvert_mod, ONLY: presnivs, ap, bp, dpres, nivsig, nivsigs, pa, & preff, aps, bps, pseudoalt, scaleheight USE temps_mod, ONLY: annee_ref, calend, day_end, day_ini, day_ref, & itau_dyn, itau_phy, start_time, year_len USE phys_cal_mod, ONLY: year_len_phys_cal_mod => year_len USE lmdz_1dutils, ONLY: fq_sat, conf_unicol, dyn1deta0, dyn1dredem, disvert0 USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_OUTPUTPHYSSCM INCLUDE "dimensions.h" INCLUDE "YOMCST.h" INCLUDE "clesphys.h" INCLUDE "dimsoil.h" INCLUDE "compar1d.h" INCLUDE "flux_arp.h" INCLUDE "date_cas.h" INCLUDE "tsoilnudge.h" INCLUDE "fcg_gcssold.h" INCLUDE "compbl.h" !===================================================================== ! DECLARATIONS !===================================================================== !--------------------------------------------------------------------- ! Arguments d' initialisations de la physique (USER DEFINE) !--------------------------------------------------------------------- INTEGER, parameter :: ngrid = 1 REAL :: zcufi = 1. REAL :: zcvfi = 1. REAL :: fnday REAL :: day, daytime REAL :: day1 REAL :: heure INTEGER :: jour INTEGER :: mois INTEGER :: an !--------------------------------------------------------------------- ! Declarations related to forcing and initial profiles !--------------------------------------------------------------------- INTEGER :: kmax = llm INTEGER llm700, nq1, nq2 INTEGER, PARAMETER :: nlev_max = 1000, nqmx = 1000 REAL timestep, frac REAL height(nlev_max), tttprof(nlev_max), qtprof(nlev_max) real uprof(nlev_max), vprof(nlev_max), e12prof(nlev_max) real ugprof(nlev_max), vgprof(nlev_max), wfls(nlev_max) real dqtdxls(nlev_max), dqtdyls(nlev_max) real dqtdtls(nlev_max), thlpcar(nlev_max) real qprof(nlev_max, nqmx) ! INTEGER :: forcing_type LOGICAL :: forcing_les = .FALSE. LOGICAL :: forcing_armcu = .FALSE. LOGICAL :: forcing_rico = .FALSE. LOGICAL :: forcing_radconv = .FALSE. LOGICAL :: forcing_toga = .FALSE. LOGICAL :: forcing_twpice = .FALSE. LOGICAL :: forcing_amma = .FALSE. LOGICAL :: forcing_dice = .FALSE. LOGICAL :: forcing_gabls4 = .FALSE. LOGICAL :: forcing_GCM2SCM = .FALSE. LOGICAL :: forcing_GCSSold = .FALSE. LOGICAL :: forcing_sandu = .FALSE. LOGICAL :: forcing_astex = .FALSE. LOGICAL :: forcing_fire = .FALSE. LOGICAL :: forcing_case = .FALSE. LOGICAL :: forcing_case2 = .FALSE. LOGICAL :: forcing_SCM = .FALSE. !flag forcings LOGICAL :: nudge_wind = .TRUE. LOGICAL :: nudge_thermo = .FALSE. LOGICAL :: cptadvw = .TRUE. !===================================================================== ! DECLARATIONS FOR EACH CASE !===================================================================== INCLUDE "1D_decl_cases.h" !--------------------------------------------------------------------- ! Declarations related to nudging !--------------------------------------------------------------------- INTEGER :: nudge_max parameter (nudge_max = 9) INTEGER :: inudge_RHT = 1 INTEGER :: inudge_UV = 2 LOGICAL :: nudge(nudge_max) REAL :: t_targ(llm) REAL :: rh_targ(llm) REAL :: u_targ(llm) REAL :: v_targ(llm) !--------------------------------------------------------------------- ! Declarations related to vertical discretization: !--------------------------------------------------------------------- REAL :: pzero = 1.e5 REAL :: play (llm), zlay (llm), sig_s(llm), plev(llm + 1) REAL :: playd(llm), zlayd(llm), ap_amma(llm + 1), bp_amma(llm + 1) !--------------------------------------------------------------------- ! Declarations related to variables !--------------------------------------------------------------------- REAL :: phi(llm) REAL :: teta(llm), tetal(llm), temp(llm), u(llm), v(llm), w(llm) REAL rot(1, llm) ! relative vorticity, in s-1 REAL :: rlat_rad(1), rlon_rad(1) REAL :: omega(llm), omega2(llm), rho(llm + 1) REAL :: ug(llm), vg(llm), fcoriolis REAL :: sfdt, cfdt REAL :: du_phys(llm), dv_phys(llm), dt_phys(llm) REAL :: w_adv(llm), z_adv(llm) REAL :: d_t_vert_adv(llm), d_u_vert_adv(llm), d_v_vert_adv(llm) REAL :: dt_cooling(llm), d_t_adv(llm), d_t_nudge(llm) REAL :: d_u_nudge(llm), d_v_nudge(llm) ! REAL :: d_u_adv(llm),d_v_adv(llm) REAL :: d_u_age(llm), d_v_age(llm) REAL :: alpha REAL :: ttt REAL, ALLOCATABLE, DIMENSION(:, :) :: q REAL, ALLOCATABLE, DIMENSION(:, :) :: dq REAL, ALLOCATABLE, DIMENSION(:, :) :: d_q_vert_adv REAL, ALLOCATABLE, DIMENSION(:, :) :: d_q_adv REAL, ALLOCATABLE, DIMENSION(:, :) :: d_q_nudge ! REAL, ALLOCATABLE, DIMENSION(:):: d_th_adv !--------------------------------------------------------------------- ! Initialization of surface variables !--------------------------------------------------------------------- REAL :: run_off_lic_0(1) REAL :: fder(1), snsrf(1, nbsrf), qsurfsrf(1, nbsrf) REAL :: tsoil(1, nsoilmx, nbsrf) ! REAL :: agesno(1,nbsrf) !--------------------------------------------------------------------- ! Call to phyredem !--------------------------------------------------------------------- LOGICAL :: ok_writedem = .TRUE. REAL :: sollw_in = 0. REAL :: solsw_in = 0. !--------------------------------------------------------------------- ! Call to physiq !--------------------------------------------------------------------- LOGICAL :: firstcall = .TRUE. LOGICAL :: lastcall = .FALSE. REAL :: phis(1) = 0.0 REAL :: dpsrf(1) !--------------------------------------------------------------------- ! Initializations of boundary conditions !--------------------------------------------------------------------- REAL, ALLOCATABLE :: phy_nat (:) ! 0=ocean libre,1=land,2=glacier,3=banquise REAL, ALLOCATABLE :: phy_alb (:) ! Albedo land only (old value condsurf_jyg=0.3) REAL, ALLOCATABLE :: phy_sst (:) ! SST (will not be used; cf read_tsurf1d.F) REAL, ALLOCATABLE :: phy_bil (:) ! Ne sert que pour les slab_ocean REAL, ALLOCATABLE :: phy_rug (:) ! Longueur rugosite utilisee sur land only REAL, ALLOCATABLE :: phy_ice (:) ! Fraction de glace REAL, ALLOCATABLE :: phy_fter(:) ! Fraction de terre REAL, ALLOCATABLE :: phy_foce(:) ! Fraction de ocean REAL, ALLOCATABLE :: phy_fsic(:) ! Fraction de glace REAL, ALLOCATABLE :: phy_flic(:) ! Fraction de glace !--------------------------------------------------------------------- ! Fichiers et d'autres variables !--------------------------------------------------------------------- INTEGER :: k, l, i, it = 1, mxcalc INTEGER :: nsrf INTEGER jcode INTEGER read_climoz INTEGER :: it_end ! iteration number of the last call !Al1,plev,play,phi,phis,presnivs, INTEGER ecrit_slab_oc !1=ecrit,-1=lit,0=no file data ecrit_slab_oc/-1/ ! if flag_inhib_forcing = 0, tendencies of forcing are added ! <> 0, tendencies of forcing are not added INTEGER :: flag_inhib_forcing = 0 PRINT*, 'VOUS ENTREZ DANS LE 1D FORMAT STANDARD' !===================================================================== ! INITIALIZATIONS !===================================================================== du_phys(:) = 0. dv_phys(:) = 0. dt_phys(:) = 0. d_t_vert_adv(:) = 0. d_u_vert_adv(:) = 0. d_v_vert_adv(:) = 0. dt_cooling(:) = 0. d_t_adv(:) = 0. d_t_nudge(:) = 0. d_u_nudge(:) = 0. d_v_nudge(:) = 0. d_u_adv(:) = 0. d_v_adv(:) = 0. d_u_age(:) = 0. d_v_age(:) = 0. ! Initialization of Common turb_forcing dtime_frcg = 0. Turb_fcg_gcssold = .FALSE. hthturb_gcssold = 0. hqturb_gcssold = 0. !--------------------------------------------------------------------- ! OPTIONS OF THE 1D SIMULATION (lmdz1d.def => unicol.def) !--------------------------------------------------------------------- CALL conf_unicol !Al1 moves this gcssold var from common fcg_gcssold to Turb_fcg_gcssold = xTurb_fcg_gcssold ! -------------------------------------------------------------------- close(1) WRITE(*, *) 'lmdz1d.def lu => unicol.def' forcing_SCM = .TRUE. year_ini_cas = 1997 ! It is possible that those parameters are run twice. ! A REVOIR : LIRE PEUT ETRE AN MOIS JOUR DIRECETEMENT CALL getin('anneeref', year_ini_cas) CALL getin('dayref', day_deb) mth_ini_cas = 1 ! pour le moment on compte depuis le debut de l'annee CALL getin('time_ini', heure_ini_cas) PRINT*, 'NATURE DE LA SURFACE ', nat_surf ! Initialization of the LOGICAL switch for nudging jcode = iflag_nudge do i = 1, nudge_max nudge(i) = mod(jcode, 10) >= 1 jcode = jcode / 10 enddo !----------------------------------------------------------------------- ! Definition of the run !----------------------------------------------------------------------- CALL conf_gcm(99, .TRUE.) !----------------------------------------------------------------------- allocate(phy_nat (year_len)) ! 0=ocean libre,1=land,2=glacier,3=banquise phy_nat(:) = 0.0 allocate(phy_alb (year_len)) ! Albedo land only (old value condsurf_jyg=0.3) allocate(phy_sst (year_len)) ! SST (will not be used; cf read_tsurf1d.F) allocate(phy_bil (year_len)) ! Ne sert que pour les slab_ocean phy_bil(:) = 1.0 allocate(phy_rug (year_len)) ! Longueur rugosite utilisee sur land only allocate(phy_ice (year_len)) ! Fraction de glace phy_ice(:) = 0.0 allocate(phy_fter(year_len)) ! Fraction de terre phy_fter(:) = 0.0 allocate(phy_foce(year_len)) ! Fraction de ocean phy_foce(:) = 0.0 allocate(phy_fsic(year_len)) ! Fraction de glace phy_fsic(:) = 0.0 allocate(phy_flic(year_len)) ! Fraction de glace phy_flic(:) = 0.0 !----------------------------------------------------------------------- ! Choix du calendrier ! ------------------- ! calend = 'earth_365d' IF (calend == 'earth_360d') THEN CALL ioconf_calendar('360_day') WRITE(*, *)'CALENDRIER CHOISI: Terrestre a 360 jours/an' ELSE IF (calend == 'earth_365d') THEN CALL ioconf_calendar('noleap') WRITE(*, *)'CALENDRIER CHOISI: Terrestre a 365 jours/an' ELSE IF (calend == 'earth_366d') THEN CALL ioconf_calendar('all_leap') WRITE(*, *)'CALENDRIER CHOISI: Terrestre bissextile' ELSE IF (calend == 'gregorian') THEN stop 'gregorian calend should not be used by normal user' CALL ioconf_calendar('gregorian') ! not to be used by normal users WRITE(*, *)'CALENDRIER CHOISI: Gregorien' else write (*, *) 'ERROR : unknown calendar ', calend stop 'calend should be 360d,earth_365d,earth_366d,gregorian' endif !----------------------------------------------------------------------- !c Date : ! La date est supposee donnee sous la forme [annee, numero du jour dans ! l annee] ; l heure est donnee dans time_ini, lu dans lmdz1d.def. ! On appelle ymds2ju pour convertir [annee, jour] en [jour Julien]. ! Le numero du jour est dans "day". L heure est traitee separement. ! La date complete est dans "daytime" (l'unite est le jour). IF (nday>0) THEN fnday = nday else fnday = -nday / float(day_step) endif print *, 'fnday=', fnday ! start_time doit etre en FRACTION DE JOUR start_time = time_ini / 24. annee_ref = anneeref mois = 1 day_ref = dayref heure = 0. itau_dyn = 0 itau_phy = 0 CALL ymds2ju(annee_ref, mois, day_ref, heure, day) day_ini = int(day) day_end = day_ini + int(fnday) ! Convert the initial date to Julian day day_ini_cas = day_deb PRINT*, 'time case', year_ini_cas, mth_ini_cas, day_ini_cas CALL ymds2ju & (year_ini_cas, mth_ini_cas, day_ini_cas, heure_ini_cas * 3600 & , day_ju_ini_cas) PRINT*, 'time case 2', day_ini_cas, day_ju_ini_cas daytime = day + heure_ini_cas / 24. ! 1st day and initial time of the simulation ! Print out the actual date of the beginning of the simulation : CALL ju2ymds(daytime, year_print, month_print, day_print, sec_print) print *, ' Time of beginning : ', & year_print, month_print, day_print, sec_print !--------------------------------------------------------------------- ! Initialization of dimensions, geometry and initial state !--------------------------------------------------------------------- ! CALL init_phys_lmdz(1,1,llm,1,(/1/)) ! job now done via iniphysiq ! but we still need to initialize dimphy module (klon,klev,etc.) here. CALL init_dimphy1D(1, llm) CALL suphel CALL init_infotrac IF (nqtot>nqmx) STOP 'Augmenter nqmx dans lmdz1d.F' allocate(q(llm, nqtot)) ; q(:, :) = 0. allocate(dq(llm, nqtot)) allocate(d_q_vert_adv(llm, nqtot)) allocate(d_q_adv(llm, nqtot)) allocate(d_q_nudge(llm, nqtot)) ! allocate(d_th_adv(llm)) q(:, :) = 0. dq(:, :) = 0. d_q_vert_adv(:, :) = 0. d_q_adv(:, :) = 0. d_q_nudge(:, :) = 0. ! No ozone climatology need be read in this pre-initialization ! (phys_state_var_init is called again in physiq) read_climoz = 0 nsw = 6 CALL phys_state_var_init(read_climoz) IF (ngrid/=klon) THEN PRINT*, 'stop in inifis' PRINT*, 'Probleme de dimensions :' PRINT*, 'ngrid = ', ngrid PRINT*, 'klon = ', klon stop endif !!!===================================================================== !!! Feedback forcing values for Gateaux differentiation (al1) !!!===================================================================== !! qsol = qsolinp qsurf = fq_sat(tsurf, psurf / 100.) beta_aridity(:, :) = beta_surf day1 = day_ini time = daytime - day ts_toga(1) = tsurf ! needed by read_tsurf1d.F rho(1) = psurf / (rd * tsurf * (1. + (rv / rd - 1.) * qsurf)) !! mpl et jyg le 22/08/2012 : !! pour que les cas a flux de surface imposes marchent IF(.NOT.ok_flux_surf.OR.max(abs(wtsurf), abs(wqsurf))>0.) THEN fsens = -wtsurf * rcpd * rho(1) flat = -wqsurf * rlvtt * rho(1) print *, 'Flux: ok_flux wtsurf wqsurf', ok_flux_surf, wtsurf, wqsurf ENDIF PRINT*, 'Flux sol ', fsens, flat ! Vertical discretization and pressure levels at half and mid levels: pa = 5e4 !! preff= 1.01325e5 preff = psurf IF (ok_old_disvert) THEN CALL disvert0(pa, preff, ap, bp, dpres, presnivs, nivsigs, nivsig) print *, 'On utilise disvert0' aps(1:llm) = 0.5 * (ap(1:llm) + ap(2:llm + 1)) bps(1:llm) = 0.5 * (bp(1:llm) + bp(2:llm + 1)) scaleheight = 8. pseudoalt(1:llm) = -scaleheight * log(presnivs(1:llm) / preff) ELSE CALL disvert() print *, 'On utilise disvert' ! Nouvelle version disvert permettant d imposer ap,bp (modif L.Guez) MPL 18092012 ! Dans ce cas, on lit ap,bp dans le fichier hybrid.txt ENDIF sig_s = presnivs / preff plev = ap + bp * psurf play = 0.5 * (plev(1:llm) + plev(2:llm + 1)) zlay = -rd * 300. * log(play / psurf) / rg ! moved after reading profiles. IF (forcing_type == 59) THEN ! pour forcing_sandu, on cherche l'indice le plus proche de 700hpa#3000m WRITE(*, *) '***********************' do l = 1, llm WRITE(*, *) 'l,play(l),presnivs(l): ', l, play(l), presnivs(l) IF (trouve_700 .AND. play(l)<=70000) THEN llm700 = l print *, 'llm700,play=', llm700, play(l) / 100. trouve_700 = .FALSE. endif enddo WRITE(*, *) '***********************' ENDIF !===================================================================== ! EVENTUALLY, READ FORCING DATA : !===================================================================== INCLUDE "1D_read_forc_cases.h" PRINT*, 'A d_t_adv ', d_t_adv(1:20)*86400 IF (forcing_GCM2SCM) THEN write (*, *) 'forcing_GCM2SCM not yet implemented' stop 'in initialization' ENDIF ! forcing_GCM2SCM !===================================================================== ! Initialisation de la physique : !===================================================================== ! Rq: conf_phys.F90 lit tous les flags de physiq.def; conf_phys appele depuis physiq.F ! day_step, iphysiq lus dans gcm.def ci-dessus ! timestep: calcule ci-dessous from rday et day_step ! ngrid=1 ! llm: defini dans .../modipsl/modeles/LMDZ4/libf/grid/dimension ! rday: defini dans suphel.F (86400.) ! day_ini: lu dans run.def (dayref) ! rlat_rad,rlon-rad: transformes en radian de rlat,rlon lus dans lmdz1d.def (en degres) ! airefi,zcufi,zcvfi initialises au debut de ce programme ! rday,ra,rg,rd,rcpd declares dans YOMCST.h et calcules dans suphel.F day_step = float(nsplit_phys)*day_step/float(iphysiq) write (*, *) 'Time step divided by nsplit_phys (=', nsplit_phys, ')' timestep = rday/day_step dtime_frcg = timestep zcufi = airefi zcvfi = airefi rlat_rad(1) = xlat*rpi/180. rlon_rad(1) = xlon*rpi/180. ! iniphysiq will CALL iniaqua who needs year_len from phys_cal_mod year_len_phys_cal_mod = year_len ! Ehouarn: iniphysiq requires arrays related to (3D) dynamics grid, ! e.g. for cell boundaries, which are meaningless in 1D; so pad these ! with '0.' when necessary CALL iniphysiq(iim, jjm, llm, & 1, comm_lmdz, & rday, day_ini, timestep, & (/rlat_rad(1), 0./), (/0./), & (/0., 0./), (/rlon_rad(1), 0./), & (/ (/airefi, 0./), (/0., 0./) /), & (/zcufi, 0., 0., 0./), & (/zcvfi, 0./), & ra, rg, rd,rcpd, 1) PRINT*, 'apres iniphysiq' ! 2 PARAMETRES QUI DEVRAIENT ETRE LUS DANS run.def MAIS NE LE SONT PAS ICI: co2_ppm = 330.0 solaire = 1370.0 ! Ecriture du startphy avant le premier appel a la physique. ! On le met juste avant pour avoir acces a tous les champs IF (ok_writedem) THEN !-------------------------------------------------------------------------- ! pbl_surface_init (called here) and pbl_surface_final (called by phyredem) ! need : qsol fder snow qsurf evap rugos agesno ftsoil !-------------------------------------------------------------------------- type_ocean = "force" run_off_lic_0(1) = restart_runoff CALL fonte_neige_init(run_off_lic_0) fder = 0. snsrf(1, :) = snowmass ! masse de neige des sous surface qsurfsrf(1, :) = qsurf ! humidite de l'air des sous surface fevap = 0. z0m(1, :) = rugos ! couverture de neige des sous surface z0h(1, :) = rugosh ! couverture de neige des sous surface agesno = xagesno tsoil(:, :, :) = tsurf !------ AMMA 2e run avec modele sol et rayonnement actif (MPL 23052012) ! tsoil(1,1,1)=299.18 ! tsoil(1,2,1)=300.08 ! tsoil(1,3,1)=301.88 ! tsoil(1,4,1)=305.48 ! tsoil(1,5,1)=308.00 ! tsoil(1,6,1)=308.00 ! tsoil(1,7,1)=308.00 ! tsoil(1,8,1)=308.00 ! tsoil(1,9,1)=308.00 ! tsoil(1,10,1)=308.00 ! tsoil(1,11,1)=308.00 !----------------------------------------------------------------------- CALL pbl_surface_init(fder, snsrf, qsurfsrf, tsoil) !------------------ prepare limit conditions for limit.nc ----------------- !-- Ocean force PRINT*, 'avant phyredem' pctsrf(1, :) = 0. IF (nat_surf==0.) THEN pctsrf(1, is_oce) = 1. pctsrf(1, is_ter) = 0. pctsrf(1, is_lic) = 0. pctsrf(1, is_sic) = 0. ELSE IF (nat_surf == 1) THEN pctsrf(1, is_oce) = 0. pctsrf(1, is_ter) = 1. pctsrf(1, is_lic) = 0. pctsrf(1, is_sic) = 0. ELSE IF (nat_surf == 2) THEN pctsrf(1, is_oce) = 0. pctsrf(1, is_ter) = 0. pctsrf(1, is_lic) = 1. pctsrf(1, is_sic) = 0. ELSE IF (nat_surf == 3) THEN pctsrf(1, is_oce) = 0. pctsrf(1, is_ter) = 0. pctsrf(1, is_lic) = 0. pctsrf(1, is_sic) = 1. end if PRINT*, 'nat_surf,pctsrf(1,is_oce),pctsrf(1,is_ter)', nat_surf & , pctsrf(1, is_oce), pctsrf(1, is_ter) zmasq = pctsrf(1, is_ter)+pctsrf(1, is_lic) zpic = zpicinp ftsol = tsurf falb_dir= albedo falb_dif = albedo rugoro = rugos t_ancien(1, :)= temp(:) q_ancien(1, :)= q(:, 1) ql_ancien = 0. qs_ancien = 0. prlw_ancien = 0. prsw_ancien = 0. prw_ancien = 0. !jyg< ! Etienne: comment those lines since now the TKE is inialized in 1D_read_forc_cases !! pbl_tke(:,:,:)=1.e-8 ! pbl_tke(:,:,:)=0. ! pbl_tke(:,2,:)=1.e-2 !>jyg rain_fall = 0. snow_fall = 0. solsw = 0. solswfdiff= 0. sollw = 0. sollwdown = rsigma*tsurf**4 radsol = 0. rnebcon= 0. ratqs = 0. clwcon = 0. zmax0 = 0. zmea = zsurf zstd= 0. zsig = 0. zgam = 0. zval = 0. zthe = 0. sig1= 0. w01 = 0. wake_deltaq = 0. wake_deltat = 0. wake_delta_pbl_TKE(:, :, :) = 0. delta_tsurf = 0. wake_fip = 0. wake_pe = 0. wake_s = 0. awake_s = 0. wake_dens = 0. awake_dens = 0. cv_gen = 0. wake_cstar = 0. ale_bl = 0. ale_bl_trig = 0. alp_bl = 0. IF (ALLOCATED(du_gwd_rando)) du_gwd_rando = 0. IF (ALLOCATED(du_gwd_front)) du_gwd_front = 0. entr_therm = 0. detr_therm = 0. f0 = 0. fm_therm = 0. u_ancien(1, :)= u(:) v_ancien(1, :)= v(:) rneb_ancien(1, :)= 0. u10m = 0. v10m = 0. ale_wake = 0. ale_bl_stat = 0. ratqs_inter_(:, :)= 0.002 !------------------------------------------------------------------------ ! Make file containing restart for the physics (startphy.nc) ! NB: List of the variables to be written by phyredem (via put_field): ! rlon,rlat,zmasq,pctsrf(:,is_ter),pctsrf(:,is_lic),pctsrf(:,is_oce) ! pctsrf(:,is_sic),ftsol(:,nsrf),tsoil(:,isoil,nsrf),qsurf(:,nsrf) ! qsol,falb_dir(:,nsrf),falb_dif(:,nsrf),evap(:,nsrf),snow(:,nsrf) ! radsol,solsw,solswfdiff,sollw, sollwdown,fder,rain_fall,snow_fall,frugs(:,nsrf) ! agesno(:,nsrf),zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro ! t_ancien,q_ancien,,frugs(:,is_oce),clwcon(:,1),rnebcon(:,1),ratqs(:,1) ! run_off_lic_0,pbl_tke(:,1:klev,nsrf), zmax0,f0,sig1,w01 ! wake_deltat,wake_deltaq,wake_s,awake_s,wake_dens,awake_dens,cv_gen,wake_cstar, ! wake_fip,wake_delta_pbl_tke(:,1:klev,nsrf) ! NB2: The content of the startphy.nc file depends on some flags defined in ! the ".def" files. However, since conf_phys is not called in lmdz1d.F90, these flags have ! to be set at some arbitratry convenient values. !------------------------------------------------------------------------ !Al1 =============== restart option ====================================== iflag_physiq = 0 CALL getin('iflag_physiq', iflag_physiq) IF (.NOT.restart) THEN iflag_pbl = 5 CALL phyredem ("startphy.nc") else ! (desallocations) PRINT*, 'callin surf final' CALL pbl_surface_final(fder, snsrf, qsurfsrf, tsoil) PRINT*, 'after surf final' CALL fonte_neige_final(run_off_lic_0) endif ok_writedem = .FALSE. PRINT*,'apres phyredem' endif ! ok_writedem !------------------------------------------------------------------------ ! Make file containing boundary conditions (limit.nc) **Al1->restartdyn*** ! -------------------------------------------------- ! NB: List of the variables to be written in limit.nc ! (by writelim.F, SUBROUTINE of 1DUTILS.h): ! phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,phy_ice, ! phy_fter,phy_foce,phy_flic,phy_fsic) !------------------------------------------------------------------------ do i = 1, year_len phy_nat(i) = nat_surf phy_alb(i) = albedo phy_sst(i) = tsurf ! read_tsurf1d will be used instead phy_rug(i) = rugos phy_fter(i) = pctsrf(1, is_ter) phy_foce(i) = pctsrf(1, is_oce) phy_fsic(i) = pctsrf(1, is_sic) phy_flic(i) = pctsrf(1, is_lic) enddo ! fabrication de limit.nc CALL writelim (1, phy_nat, phy_alb, phy_sst, phy_bil,phy_rug, & phy_ice, phy_fter, phy_foce, phy_flic,phy_fsic) CALL phys_state_var_end !Al1 IF (restart) THEN PRINT*, 'CALL to restart dyn 1d' Call dyn1deta0("start1dyn.nc", plev, play, phi, phis,presnivs, & u, v, temp, q,omega2) PRINT*, 'fnday,annee_ref,day_ref,day_ini', & fnday, annee_ref,day_ref, day_ini !** CALL ymds2ju(annee_ref,mois,day_ini,heure,day) day = day_ini day_end = day_ini + nday daytime = day + time_ini/24. ! 1st day and initial time of the simulation ! Print out the actual date of the beginning of the simulation : CALL ju2ymds(daytime, an, mois, jour, heure) print *, ' Time of beginning : y m d h', an, mois,jour, heure/3600. day = int(daytime) time = daytime-day PRINT*,'****** intialised fields from restart1dyn *******' PRINT*, 'plev,play,phi,phis,presnivs,u,v,temp,q,omega2' PRINT*,'temp(1),q(1,1),u(1),v(1),plev(1),phis :' PRINT*, temp(1), q(1, 1), u(1), v(1), plev(1), phis(1) ! raz for safety do l = 1, llm d_q_vert_adv(l, 1) = 0. enddo endif !====================== end restart ================================= IF (ecrit_slab_oc==1) THEN open(97, file = 'div_slab.dat', STATUS = 'UNKNOWN') elseif (ecrit_slab_oc==0) THEN open(97, file = 'div_slab.dat', STATUS = 'OLD') endif !===================================================================== IF (CPPKEY_OUTPUTPHYSSCM) THEN CALL iophys_ini(timestep) END IF !===================================================================== ! START OF THE TEMPORAL LOOP : !===================================================================== it_end = nint(fnday*day_step) do while(it<=it_end) IF (prt_level>=1) THEN PRINT*, 'XXXXXXXXXXXXXXXXXXX ITAP,day,time=', & it, day, time, it_end, day_step PRINT*,'PAS DE TEMPS ', timestep endif IF (it==it_end) lastcall = .True. !--------------------------------------------------------------------- ! Interpolation of forcings in time and onto model levels !--------------------------------------------------------------------- INCLUDE "1D_interp_cases.h" !--------------------------------------------------------------------- ! Geopotential : !--------------------------------------------------------------------- phis(1)= zsurf*RG ! phi(1)=phis(1)+RD*temp(1)*(plev(1)-play(1))/(.5*(plev(1)+play(1))) ! Calculate geopotential from the ground surface since phi and phis are added in physiq_mod phi(1)= RD*temp(1)*(plev(1)-play(1))/(.5*(plev(1)+play(1))) do l = 1, llm-1 phi(l+1)= phi(l)+RD*(temp(l)+temp(l+1))* & (play(l)-play(l+1))/(play(l)+play(l+1)) enddo !--------------------------------------------------------------------- ! Vertical advection !--------------------------------------------------------------------- IF (forc_w+forc_omega > 0) THEN IF (forc_w == 1) THEN w_adv = w_mod_cas z_adv = phi/RG ELSE w_adv = omega z_adv =play ENDIF teta = temp*(pzero/play)**rkappa do l = 2, llm-1 ! vertical tendencies computed as d X / d t = -W d X / d z d_u_vert_adv(l)= -w_adv(l)*(u(l+1)-u(l-1))/(z_adv(l+1)-z_adv(l-1)) d_v_vert_adv(l)= -w_adv(l)*(v(l+1)-v(l-1))/(z_adv(l+1)-z_adv(l-1)) ! d theta / dt = -W d theta / d z, transformed into d temp / d t dividing by (pzero/play(l))**rkappa d_t_vert_adv(l)= -w_adv(l)*(teta(l+1)-teta(l-1))/(z_adv(l+1)-z_adv(l-1)) / (pzero/play(l))**rkappa d_q_vert_adv(l, 1)= -w_adv(l)*(q(l+1, 1)-q(l-1, 1))/(z_adv(l+1)-z_adv(l-1)) enddo d_u_adv(:)= d_u_adv(:)+d_u_vert_adv(:) d_v_adv(:)= d_v_adv(:)+d_v_vert_adv(:) d_t_adv(:)= d_t_adv(:)+d_t_vert_adv(:) d_q_adv(:, 1)= d_q_adv(:, 1)+d_q_vert_adv(:, 1) ENDIF !--------------------------------------------------------------------- ! Listing output for debug prt_level>=1 !--------------------------------------------------------------------- IF (prt_level>=1) THEN print *, ' avant physiq : -------- day time ', day, time WRITE(*, *) 'firstcall,lastcall,phis', & firstcall, lastcall, phis end if IF (prt_level>=5) THEN WRITE(*, '(a10,2a4,4a13)') 'BEFOR1 IT=', 'it', 'l', & 'presniv', 'plev','play', 'phi' WRITE(*, '(a10,2i4,4f13.2)') ('BEFOR1 IT= ', it, l, & presnivs(l), plev(l), play(l), phi(l), l = 1, llm) WRITE(*, '(a11,2a4,a11,6a8)') 'BEFOR2', 'it', 'l', & 'presniv', 'u','v', 'temp', 'q1', 'q2', 'omega2' WRITE(*, '(a11,2i4,f11.2,5f8.2,e10.2)') ('BEFOR2 IT= ', it, l, & presnivs(l), u(l), v(l), temp(l), q(l, 1), q(l, 2), omega2(l), l = 1, llm) endif !--------------------------------------------------------------------- ! Call physiq : !--------------------------------------------------------------------- CALL physiq(ngrid, llm, & firstcall, lastcall, timestep, & plev, play, phi, phis, presnivs, & u, v, rot, temp, q,omega2, & du_phys, dv_phys, dt_phys, dq,dpsrf) firstcall = .FALSE. !--------------------------------------------------------------------- ! Listing output for debug !--------------------------------------------------------------------- IF (prt_level>=5) THEN WRITE(*, '(a11,2a4,4a13)') 'AFTER1 IT=', 'it', 'l', & 'presniv', 'plev','play', 'phi' WRITE(*, '(a11,2i4,4f13.2)') ('AFTER1 it= ', it, l, & presnivs(l), plev(l), play(l), phi(l), l = 1, llm) WRITE(*, '(a11,2a4,a11,6a8)') 'AFTER2', 'it', 'l', & 'presniv', 'u','v', 'temp', 'q1', 'q2', 'omega2' WRITE(*, '(a11,2i4,f11.2,5f8.2,e10.2)') ('AFTER2 it= ', it, l, & presnivs(l), u(l), v(l), temp(l), q(l, 1), q(l, 2), omega2(l), l = 1, llm) WRITE(*, '(a11,2a4,a11,5a8)') 'AFTER3', 'it', 'l', & 'presniv', 'du_phys','dv_phys', 'dt_phys', 'dq1', 'dq2' WRITE(*, '(a11,2i4,f11.2,5f8.2)') ('AFTER3 it= ', it, l, & presnivs(l), 86400*du_phys(l), 86400*dv_phys(l), & 86400*dt_phys(l), 86400*dq(l, 1), dq(l, 2), l = 1, llm) WRITE(*, *) 'dpsrf', dpsrf endif !--------------------------------------------------------------------- ! Add physical tendencies : !--------------------------------------------------------------------- fcoriolis= 2.*sin(rpi*xlat/180.)*romega IF (prt_level >= 5) PRINT*, 'fcoriolis, xlat,mxcalc ', & fcoriolis, xlat, mxcalc !--------------------------------------------------------------------- ! Geostrophic forcing !--------------------------------------------------------------------- IF (forc_geo == 0) THEN d_u_age(1:mxcalc)= 0. d_v_age(1:mxcalc)= 0. ELSE sfdt = sin(0.5*fcoriolis*timestep) cfdt = cos(0.5*fcoriolis*timestep) d_u_age(1:mxcalc)= -2.*sfdt/timestep* & (sfdt*(u(1:mxcalc)-ug(1:mxcalc)) - & cfdt*(v(1:mxcalc)-vg(1:mxcalc))) !! : fcoriolis*(v(1:mxcalc)-vg(1:mxcalc)) d_v_age(1:mxcalc)= -2.*sfdt/timestep* & (cfdt*(u(1:mxcalc)-ug(1:mxcalc)) + & sfdt*(v(1:mxcalc)-vg(1:mxcalc))) !! : -fcoriolis*(u(1:mxcalc)-ug(1:mxcalc)) ENDIF !--------------------------------------------------------------------- ! Nudging ! EV: rewrite the section to account for a time- and height-varying ! nudging !--------------------------------------------------------------------- d_t_nudge(:) = 0. d_u_nudge(:) = 0. d_v_nudge(:) = 0. d_q_nudge(:, :) = 0. DO l = 1, llm IF (nudging_u < 0) THEN d_u_nudge(l)=(u_nudg_mod_cas(l)-u(l))*invtau_u_nudg_mod_cas(l) ELSE IF (play(l) < p_nudging_u .AND. nint(nudging_u) /= 0) & d_u_nudge(l)=(u_nudg_mod_cas(l)-u(l))/nudging_u ENDIF IF (nudging_v < 0) THEN d_v_nudge(l)=(v_nudg_mod_cas(l)-v(l))*invtau_v_nudg_mod_cas(l) ELSE IF (play(l) < p_nudging_v .AND. nint(nudging_v) /= 0) & d_v_nudge(l)=(v_nudg_mod_cas(l)-v(l))/nudging_v ENDIF IF (nudging_t < 0) THEN d_t_nudge(l)=(temp_nudg_mod_cas(l)-temp(l))*invtau_temp_nudg_mod_cas(l) ELSE IF (play(l) < p_nudging_t .AND. nint(nudging_t) /= 0) & d_t_nudge(l)=(temp_nudg_mod_cas(l)-temp(l))/nudging_t ENDIF IF (nudging_qv < 0) THEN d_q_nudge(l, 1)=(qv_nudg_mod_cas(l)-q(l, 1))*invtau_qv_nudg_mod_cas(l) ELSE IF (play(l) < p_nudging_qv .AND. nint(nudging_qv) /= 0) & d_q_nudge(l, 1)=(qv_nudg_mod_cas(l)-q(l, 1))/nudging_qv ENDIF ENDDO !--------------------------------------------------------------------- ! Optional outputs !--------------------------------------------------------------------- IF (CPPKEY_OUTPUTPHYSSCM) THEN CALL iophys_ecrit('w_adv', klev, 'w_adv', 'K/day', w_adv) CALL iophys_ecrit('z_adv', klev, 'z_adv', 'K/day', z_adv) CALL iophys_ecrit('dtadv', klev, 'dtadv', 'K/day', 86400*d_t_adv) CALL iophys_ecrit('dtdyn', klev, 'dtdyn', 'K/day', 86400*d_t_vert_adv) CALL iophys_ecrit('qv', klev, 'qv', 'g/kg', 1000*q(:, 1)) CALL iophys_ecrit('qvnud', klev, 'qvnud', 'g/kg', 1000*u_nudg_mod_cas) CALL iophys_ecrit('u', klev, 'u', 'm/s', u) CALL iophys_ecrit('unud', klev, 'unud', 'm/s', u_nudg_mod_cas) CALL iophys_ecrit('v', klev, 'v', 'm/s', v) CALL iophys_ecrit('vnud', klev, 'vnud', 'm/s', v_nudg_mod_cas) CALL iophys_ecrit('temp', klev, 'temp', 'K', temp) CALL iophys_ecrit('tempnud', klev, 'temp_nudg_mod_cas', 'K', temp_nudg_mod_cas) CALL iophys_ecrit('dtnud', klev, 'dtnud', 'K/day', 86400*d_t_nudge) CALL iophys_ecrit('dqnud', klev, 'dqnud', 'K/day', 1000*86400*d_q_nudge(:, 1)) END IF IF (flag_inhib_forcing == 0) then ! if tendency of forcings should be added u(1:mxcalc)= u(1:mxcalc) + timestep*(& du_phys(1:mxcalc) & +d_u_age(1:mxcalc)+d_u_adv(1:mxcalc) & +d_u_nudge(1:mxcalc)) v(1:mxcalc)= v(1:mxcalc) + timestep*(& dv_phys(1:mxcalc) & +d_v_age(1:mxcalc)+d_v_adv(1:mxcalc) & +d_v_nudge(1:mxcalc)) q(1:mxcalc, :)= q(1:mxcalc, :)+timestep*(& dq(1:mxcalc, :) & +d_q_adv(1:mxcalc, :) & +d_q_nudge(1:mxcalc, :)) IF (prt_level>=3) THEN print *, & 'physiq-> temp(1),dt_phys(1),d_t_adv(1),dt_cooling(1) ', & temp(1), dt_phys(1), d_t_adv(1), dt_cooling(1) PRINT*, 'dv_phys=', dv_phys PRINT* , 'd_v_age=', d_v_age PRINT*, 'd_v_adv=',d_v_adv PRINT*, 'd_v_nudge=', d_v_nudge PRINT*, v PRINT*, vg endif temp(1:mxcalc)= temp(1:mxcalc)+timestep*(& dt_phys(1:mxcalc) & +d_t_adv(1:mxcalc) & +d_t_nudge(1:mxcalc) & +dt_cooling(1:mxcalc)) ! Taux de chauffage ou refroid. !======================================================================= !! CONSERVE EN ATTENDANT QUE LE CAS EN QUESTION FONCTIONNE EN STD !! !======================================================================= teta = temp*(pzero/play)**rkappa !--------------------------------------------------------------------- ! Nudge soil temperature if requested !--------------------------------------------------------------------- IF (nudge_tsoil .AND. .NOT. lastcall) THEN ftsoil(1, isoil_nudge, :) = ftsoil(1, isoil_nudge, :) & -timestep/tau_soil_nudge*(ftsoil(1, isoil_nudge, :)-Tsoil_nudge) ENDIF !--------------------------------------------------------------------- ! Add large-scale tendencies (advection, etc) : !--------------------------------------------------------------------- !cc nrlmd !cc tmpvar=teta !cc CALL advect_vert(llm,omega,timestep,tmpvar,plev) !cc !cc teta(1:mxcalc)=tmpvar(1:mxcalc) !cc tmpvar(:)=q(:,1) !cc CALL advect_vert(llm,omega,timestep,tmpvar,plev) !cc q(1:mxcalc,1)=tmpvar(1:mxcalc) !cc tmpvar(:)=q(:,2) !cc CALL advect_vert(llm,omega,timestep,tmpvar,plev) !cc q(1:mxcalc,2)=tmpvar(1:mxcalc) END IF ! end if tendency of tendency should be added !--------------------------------------------------------------------- ! Air temperature : !--------------------------------------------------------------------- IF (lastcall) THEN PRINT*, 'Pas de temps final ', it CALL ju2ymds(daytime, an, mois, jour, heure) PRINT*, 'a la date : a m j h', an, mois, jour, heure/3600. endif ! incremente day time daytime = daytime+1./day_step day = int(daytime+0.1/day_step) ! time = max(daytime-day,0.0) !Al1&jyg: correction de bug !cc time = real(mod(it,day_step))/day_step time = time_ini/24.+real(mod(it, day_step))/day_step it = it+1 enddo IF (ecrit_slab_oc/=-1) close(97) !Al1 Call to 1D equivalent of dynredem (an,mois,jour,heure ?) ! --------------------------------------------------------------------------- CALL dyn1dredem("restart1dyn.nc", & plev, play, phi, phis,presnivs, & u, v, temp, q,omega2) CALL abort_gcm ('lmdz1d ', 'The End ', 0) END SUBROUTINE scm END MODULE lmdz_scm