[3541] | 1 | SUBROUTINE scm |
---|
| 2 | |
---|
| 3 | USE ioipsl, only: ju2ymds, ymds2ju, ioconf_calendar,getin |
---|
| 4 | USE phys_state_var_mod, ONLY : phys_state_var_init, phys_state_var_end, & |
---|
| 5 | clwcon, detr_therm, & |
---|
| 6 | qsol, fevap, z0m, z0h, agesno, & |
---|
| 7 | du_gwd_rando, du_gwd_front, entr_therm, f0, fm_therm, & |
---|
| 8 | falb_dir, falb_dif, & |
---|
| 9 | ftsol, pbl_tke, pctsrf, radsol, rain_fall, snow_fall, ratqs, & |
---|
| 10 | rnebcon, rugoro, sig1, w01, solaire_etat0, sollw, sollwdown, & |
---|
| 11 | solsw, t_ancien, q_ancien, u_ancien, v_ancien, wake_cstar, & |
---|
| 12 | wake_delta_pbl_TKE, delta_tsurf, wake_fip, wake_pe, & |
---|
| 13 | wake_deltaq, wake_deltat, wake_s, wake_dens, & |
---|
| 14 | zgam, zmax0, zmea, zpic, zsig, & |
---|
| 15 | zstd, zthe, zval, ale_bl, ale_bl_trig, alp_bl, ql_ancien, qs_ancien, & |
---|
[3592] | 16 | prlw_ancien, prsw_ancien, prw_ancien, & |
---|
| 17 | u10m,v10m,ale_wake,ale_bl_stat |
---|
| 18 | |
---|
[3541] | 19 | |
---|
| 20 | USE dimphy |
---|
| 21 | USE surface_data, only : type_ocean,ok_veget |
---|
| 22 | USE pbl_surface_mod, only : ftsoil, pbl_surface_init, & |
---|
| 23 | pbl_surface_final |
---|
| 24 | USE fonte_neige_mod, only : fonte_neige_init, fonte_neige_final |
---|
| 25 | |
---|
| 26 | USE infotrac ! new |
---|
| 27 | USE control_mod |
---|
| 28 | USE indice_sol_mod |
---|
| 29 | USE phyaqua_mod |
---|
| 30 | ! USE mod_1D_cases_read |
---|
| 31 | USE mod_1D_cases_read_std |
---|
| 32 | !USE mod_1D_amma_read |
---|
| 33 | USE print_control_mod, ONLY: lunout, prt_level |
---|
| 34 | USE iniphysiq_mod, ONLY: iniphysiq |
---|
| 35 | USE mod_const_mpi, ONLY: comm_lmdz |
---|
| 36 | USE physiq_mod, ONLY: physiq |
---|
| 37 | USE comvert_mod, ONLY: presnivs, ap, bp, dpres,nivsig, nivsigs, pa, & |
---|
| 38 | preff, aps, bps, pseudoalt, scaleheight |
---|
| 39 | USE temps_mod, ONLY: annee_ref, calend, day_end, day_ini, day_ref, & |
---|
| 40 | itau_dyn, itau_phy, start_time, year_len |
---|
| 41 | USE phys_cal_mod, ONLY : year_len_phys_cal_mod => year_len |
---|
| 42 | |
---|
| 43 | implicit none |
---|
| 44 | #include "dimensions.h" |
---|
| 45 | #include "YOMCST.h" |
---|
| 46 | !!#include "control.h" |
---|
| 47 | #include "clesphys.h" |
---|
| 48 | #include "dimsoil.h" |
---|
| 49 | !#include "indicesol.h" |
---|
| 50 | |
---|
[3686] | 51 | #include "compar1d.h" |
---|
[3541] | 52 | #include "flux_arp.h" |
---|
| 53 | #include "date_cas.h" |
---|
| 54 | #include "tsoilnudge.h" |
---|
| 55 | #include "fcg_gcssold.h" |
---|
| 56 | #include "compbl.h" |
---|
| 57 | |
---|
| 58 | !===================================================================== |
---|
| 59 | ! DECLARATIONS |
---|
| 60 | !===================================================================== |
---|
| 61 | |
---|
[3595] | 62 | #undef OUTPUT_PHYS_SCM |
---|
| 63 | |
---|
[3541] | 64 | !--------------------------------------------------------------------- |
---|
| 65 | ! Externals |
---|
| 66 | !--------------------------------------------------------------------- |
---|
| 67 | external fq_sat |
---|
| 68 | real fq_sat |
---|
| 69 | |
---|
| 70 | !--------------------------------------------------------------------- |
---|
| 71 | ! Arguments d' initialisations de la physique (USER DEFINE) |
---|
| 72 | !--------------------------------------------------------------------- |
---|
| 73 | |
---|
| 74 | integer, parameter :: ngrid=1 |
---|
| 75 | real :: zcufi = 1. |
---|
| 76 | real :: zcvfi = 1. |
---|
| 77 | |
---|
| 78 | !- real :: nat_surf |
---|
| 79 | !- logical :: ok_flux_surf |
---|
| 80 | !- real :: fsens |
---|
| 81 | !- real :: flat |
---|
| 82 | !- real :: tsurf |
---|
| 83 | !- real :: rugos |
---|
| 84 | !- real :: qsol(1:2) |
---|
| 85 | !- real :: qsurf |
---|
| 86 | !- real :: psurf |
---|
| 87 | !- real :: zsurf |
---|
| 88 | !- real :: albedo |
---|
| 89 | !- |
---|
| 90 | !- real :: time = 0. |
---|
| 91 | !- real :: time_ini |
---|
| 92 | !- real :: xlat |
---|
| 93 | !- real :: xlon |
---|
| 94 | !- real :: wtsurf |
---|
| 95 | !- real :: wqsurf |
---|
| 96 | !- real :: restart_runoff |
---|
| 97 | !- real :: xagesno |
---|
| 98 | !- real :: qsolinp |
---|
| 99 | !- real :: zpicinp |
---|
| 100 | !- |
---|
| 101 | real :: fnday |
---|
| 102 | real :: day, daytime |
---|
| 103 | real :: day1 |
---|
| 104 | real :: heure |
---|
| 105 | integer :: jour |
---|
| 106 | integer :: mois |
---|
| 107 | integer :: an |
---|
| 108 | |
---|
| 109 | !--------------------------------------------------------------------- |
---|
| 110 | ! Declarations related to forcing and initial profiles |
---|
| 111 | !--------------------------------------------------------------------- |
---|
| 112 | |
---|
| 113 | integer :: kmax = llm |
---|
| 114 | integer llm700,nq1,nq2 |
---|
| 115 | INTEGER, PARAMETER :: nlev_max=1000, nqmx=1000 |
---|
| 116 | real timestep, frac |
---|
| 117 | real height(nlev_max),tttprof(nlev_max),qtprof(nlev_max) |
---|
| 118 | real uprof(nlev_max),vprof(nlev_max),e12prof(nlev_max) |
---|
| 119 | real ugprof(nlev_max),vgprof(nlev_max),wfls(nlev_max) |
---|
| 120 | real dqtdxls(nlev_max),dqtdyls(nlev_max) |
---|
| 121 | real dqtdtls(nlev_max),thlpcar(nlev_max) |
---|
| 122 | real qprof(nlev_max,nqmx) |
---|
| 123 | |
---|
| 124 | ! integer :: forcing_type |
---|
| 125 | logical :: forcing_les = .false. |
---|
| 126 | logical :: forcing_armcu = .false. |
---|
| 127 | logical :: forcing_rico = .false. |
---|
| 128 | logical :: forcing_radconv = .false. |
---|
| 129 | logical :: forcing_toga = .false. |
---|
| 130 | logical :: forcing_twpice = .false. |
---|
| 131 | logical :: forcing_amma = .false. |
---|
| 132 | logical :: forcing_dice = .false. |
---|
| 133 | logical :: forcing_gabls4 = .false. |
---|
| 134 | |
---|
| 135 | logical :: forcing_GCM2SCM = .false. |
---|
| 136 | logical :: forcing_GCSSold = .false. |
---|
| 137 | logical :: forcing_sandu = .false. |
---|
| 138 | logical :: forcing_astex = .false. |
---|
| 139 | logical :: forcing_fire = .false. |
---|
| 140 | logical :: forcing_case = .false. |
---|
| 141 | logical :: forcing_case2 = .false. |
---|
| 142 | logical :: forcing_SCM = .false. |
---|
| 143 | integer :: type_ts_forcing ! 0 = SST constant; 1 = SST read from a file |
---|
| 144 | ! (cf read_tsurf1d.F) |
---|
| 145 | |
---|
| 146 | !flag forcings |
---|
| 147 | logical :: nudge_wind=.true. |
---|
| 148 | logical :: nudge_thermo=.false. |
---|
| 149 | logical :: cptadvw=.true. |
---|
| 150 | !===================================================================== |
---|
| 151 | ! DECLARATIONS FOR EACH CASE |
---|
| 152 | !===================================================================== |
---|
| 153 | ! |
---|
| 154 | #include "1D_decl_cases.h" |
---|
| 155 | ! |
---|
| 156 | !--------------------------------------------------------------------- |
---|
| 157 | ! Declarations related to nudging |
---|
| 158 | !--------------------------------------------------------------------- |
---|
| 159 | integer :: nudge_max |
---|
| 160 | parameter (nudge_max=9) |
---|
| 161 | integer :: inudge_RHT=1 |
---|
| 162 | integer :: inudge_UV=2 |
---|
| 163 | logical :: nudge(nudge_max) |
---|
| 164 | real :: t_targ(llm) |
---|
| 165 | real :: rh_targ(llm) |
---|
| 166 | real :: u_targ(llm) |
---|
| 167 | real :: v_targ(llm) |
---|
| 168 | ! |
---|
| 169 | !--------------------------------------------------------------------- |
---|
| 170 | ! Declarations related to vertical discretization: |
---|
| 171 | !--------------------------------------------------------------------- |
---|
| 172 | real :: pzero=1.e5 |
---|
| 173 | real :: play (llm),zlay (llm),sig_s(llm),plev(llm+1) |
---|
| 174 | real :: playd(llm),zlayd(llm),ap_amma(llm+1),bp_amma(llm+1) |
---|
| 175 | |
---|
| 176 | !--------------------------------------------------------------------- |
---|
| 177 | ! Declarations related to variables |
---|
| 178 | !--------------------------------------------------------------------- |
---|
| 179 | |
---|
| 180 | real :: phi(llm) |
---|
| 181 | real :: teta(llm),tetal(llm),temp(llm),u(llm),v(llm),w(llm) |
---|
| 182 | REAL rot(1, llm) ! relative vorticity, in s-1 |
---|
| 183 | real :: rlat_rad(1),rlon_rad(1) |
---|
| 184 | real :: omega(llm),omega2(llm),rho(llm+1) |
---|
| 185 | real :: ug(llm),vg(llm),fcoriolis |
---|
| 186 | real :: sfdt, cfdt |
---|
| 187 | real :: du_phys(llm),dv_phys(llm),dt_phys(llm) |
---|
[3595] | 188 | real :: w_adv(llm),z_adv(llm) |
---|
| 189 | real :: d_t_vert_adv(llm),d_u_vert_adv(llm),d_v_vert_adv(llm) |
---|
[3541] | 190 | real :: dt_cooling(llm),d_t_adv(llm),d_t_nudge(llm) |
---|
| 191 | real :: d_u_nudge(llm),d_v_nudge(llm) |
---|
[3693] | 192 | ! real :: d_u_adv(llm),d_v_adv(llm) |
---|
| 193 | real :: d_u_age(llm),d_v_age(llm) |
---|
[3541] | 194 | real :: alpha |
---|
| 195 | real :: ttt |
---|
| 196 | |
---|
| 197 | REAL, ALLOCATABLE, DIMENSION(:,:):: q |
---|
| 198 | REAL, ALLOCATABLE, DIMENSION(:,:):: dq |
---|
[3595] | 199 | REAL, ALLOCATABLE, DIMENSION(:,:):: d_q_vert_adv |
---|
[3541] | 200 | REAL, ALLOCATABLE, DIMENSION(:,:):: d_q_adv |
---|
| 201 | REAL, ALLOCATABLE, DIMENSION(:,:):: d_q_nudge |
---|
| 202 | ! REAL, ALLOCATABLE, DIMENSION(:):: d_th_adv |
---|
| 203 | |
---|
| 204 | !--------------------------------------------------------------------- |
---|
| 205 | ! Initialization of surface variables |
---|
| 206 | !--------------------------------------------------------------------- |
---|
| 207 | real :: run_off_lic_0(1) |
---|
| 208 | real :: fder(1),snsrf(1,nbsrf),qsurfsrf(1,nbsrf) |
---|
| 209 | real :: tsoil(1,nsoilmx,nbsrf) |
---|
| 210 | ! real :: agesno(1,nbsrf) |
---|
| 211 | |
---|
| 212 | !--------------------------------------------------------------------- |
---|
| 213 | ! Call to phyredem |
---|
| 214 | !--------------------------------------------------------------------- |
---|
| 215 | logical :: ok_writedem =.true. |
---|
| 216 | real :: sollw_in = 0. |
---|
| 217 | real :: solsw_in = 0. |
---|
| 218 | |
---|
| 219 | !--------------------------------------------------------------------- |
---|
| 220 | ! Call to physiq |
---|
| 221 | !--------------------------------------------------------------------- |
---|
| 222 | logical :: firstcall=.true. |
---|
| 223 | logical :: lastcall=.false. |
---|
| 224 | real :: phis(1) = 0.0 |
---|
| 225 | real :: dpsrf(1) |
---|
| 226 | |
---|
| 227 | !--------------------------------------------------------------------- |
---|
| 228 | ! Initializations of boundary conditions |
---|
| 229 | !--------------------------------------------------------------------- |
---|
| 230 | real, allocatable :: phy_nat (:) ! 0=ocean libre,1=land,2=glacier,3=banquise |
---|
| 231 | real, allocatable :: phy_alb (:) ! Albedo land only (old value condsurf_jyg=0.3) |
---|
| 232 | real, allocatable :: phy_sst (:) ! SST (will not be used; cf read_tsurf1d.F) |
---|
| 233 | real, allocatable :: phy_bil (:) ! Ne sert que pour les slab_ocean |
---|
| 234 | real, allocatable :: phy_rug (:) ! Longueur rugosite utilisee sur land only |
---|
| 235 | real, allocatable :: phy_ice (:) ! Fraction de glace |
---|
| 236 | real, allocatable :: phy_fter(:) ! Fraction de terre |
---|
| 237 | real, allocatable :: phy_foce(:) ! Fraction de ocean |
---|
| 238 | real, allocatable :: phy_fsic(:) ! Fraction de glace |
---|
| 239 | real, allocatable :: phy_flic(:) ! Fraction de glace |
---|
| 240 | |
---|
| 241 | !--------------------------------------------------------------------- |
---|
| 242 | ! Fichiers et d'autres variables |
---|
| 243 | !--------------------------------------------------------------------- |
---|
| 244 | integer :: k,l,i,it=1,mxcalc |
---|
| 245 | integer :: nsrf |
---|
| 246 | integer jcode |
---|
| 247 | INTEGER read_climoz |
---|
| 248 | ! |
---|
| 249 | integer :: it_end ! iteration number of the last call |
---|
| 250 | !Al1 |
---|
| 251 | integer ecrit_slab_oc !1=ecrit,-1=lit,0=no file |
---|
| 252 | data ecrit_slab_oc/-1/ |
---|
| 253 | ! |
---|
| 254 | ! if flag_inhib_forcing = 0, tendencies of forcing are added |
---|
| 255 | ! <> 0, tendencies of forcing are not added |
---|
| 256 | INTEGER :: flag_inhib_forcing = 0 |
---|
| 257 | |
---|
| 258 | |
---|
| 259 | print*,'VOUS ENTREZ DANS LE 1D FORMAT STANDARD' |
---|
| 260 | |
---|
| 261 | !===================================================================== |
---|
| 262 | ! INITIALIZATIONS |
---|
| 263 | !===================================================================== |
---|
| 264 | du_phys(:)=0. |
---|
| 265 | dv_phys(:)=0. |
---|
| 266 | dt_phys(:)=0. |
---|
[3595] | 267 | d_t_vert_adv(:)=0. |
---|
| 268 | d_u_vert_adv(:)=0. |
---|
| 269 | d_v_vert_adv(:)=0. |
---|
[3541] | 270 | dt_cooling(:)=0. |
---|
| 271 | d_t_adv(:)=0. |
---|
| 272 | d_t_nudge(:)=0. |
---|
| 273 | d_u_nudge(:)=0. |
---|
| 274 | d_v_nudge(:)=0. |
---|
[3693] | 275 | d_u_adv(:)=0. |
---|
| 276 | d_v_adv(:)=0. |
---|
| 277 | d_u_age(:)=0. |
---|
| 278 | d_v_age(:)=0. |
---|
[3541] | 279 | |
---|
| 280 | ! Initialization of Common turb_forcing |
---|
| 281 | dtime_frcg = 0. |
---|
| 282 | Turb_fcg_gcssold=.false. |
---|
| 283 | hthturb_gcssold = 0. |
---|
| 284 | hqturb_gcssold = 0. |
---|
| 285 | |
---|
| 286 | |
---|
| 287 | |
---|
| 288 | |
---|
| 289 | !--------------------------------------------------------------------- |
---|
| 290 | ! OPTIONS OF THE 1D SIMULATION (lmdz1d.def => unicol.def) |
---|
| 291 | !--------------------------------------------------------------------- |
---|
| 292 | !Al1 |
---|
[3686] | 293 | call conf_unicol |
---|
[3541] | 294 | !Al1 moves this gcssold var from common fcg_gcssold to |
---|
| 295 | Turb_fcg_gcssold = xTurb_fcg_gcssold |
---|
| 296 | ! -------------------------------------------------------------------- |
---|
| 297 | close(1) |
---|
| 298 | !Al1 |
---|
| 299 | write(*,*) 'lmdz1d.def lu => unicol.def' |
---|
| 300 | |
---|
| 301 | forcing_SCM = .true. |
---|
| 302 | year_ini_cas=1997 |
---|
| 303 | ! It is possible that those parameters are run twice. |
---|
| 304 | |
---|
| 305 | ! A REVOIR : LIRE PEUT ETRE AN MOIS JOUR DIRECETEMENT |
---|
| 306 | call getin('anneeref',year_ini_cas) |
---|
| 307 | call getin('dayref',day_deb) |
---|
| 308 | mth_ini_cas=1 ! pour le moment on compte depuis le debut de l'annee |
---|
| 309 | call getin('time_ini',heure_ini_cas) |
---|
| 310 | |
---|
| 311 | type_ts_forcing = 0 |
---|
| 312 | IF (nat_surf==0) type_ts_forcing=1 ! SST forcee sur OCEAN |
---|
| 313 | print*,'NATURE DE LA SURFACE ',nat_surf |
---|
| 314 | ! |
---|
| 315 | ! Initialization of the logical switch for nudging |
---|
| 316 | jcode = iflag_nudge |
---|
| 317 | do i = 1,nudge_max |
---|
| 318 | nudge(i) = mod(jcode,10) .ge. 1 |
---|
| 319 | jcode = jcode/10 |
---|
| 320 | enddo |
---|
| 321 | !--------------------------------------------------------------------- |
---|
| 322 | ! Definition of the run |
---|
| 323 | !--------------------------------------------------------------------- |
---|
| 324 | |
---|
| 325 | call conf_gcm( 99, .TRUE. ) |
---|
| 326 | |
---|
| 327 | !----------------------------------------------------------------------- |
---|
| 328 | allocate( phy_nat (year_len)) ! 0=ocean libre,1=land,2=glacier,3=banquise |
---|
| 329 | phy_nat(:)=0.0 |
---|
| 330 | allocate( phy_alb (year_len)) ! Albedo land only (old value condsurf_jyg=0.3) |
---|
| 331 | allocate( phy_sst (year_len)) ! SST (will not be used; cf read_tsurf1d.F) |
---|
| 332 | allocate( phy_bil (year_len)) ! Ne sert que pour les slab_ocean |
---|
| 333 | phy_bil(:)=1.0 |
---|
| 334 | allocate( phy_rug (year_len)) ! Longueur rugosite utilisee sur land only |
---|
| 335 | allocate( phy_ice (year_len)) ! Fraction de glace |
---|
| 336 | phy_ice(:)=0.0 |
---|
| 337 | allocate( phy_fter(year_len)) ! Fraction de terre |
---|
| 338 | phy_fter(:)=0.0 |
---|
| 339 | allocate( phy_foce(year_len)) ! Fraction de ocean |
---|
| 340 | phy_foce(:)=0.0 |
---|
| 341 | allocate( phy_fsic(year_len)) ! Fraction de glace |
---|
| 342 | phy_fsic(:)=0.0 |
---|
| 343 | allocate( phy_flic(year_len)) ! Fraction de glace |
---|
| 344 | phy_flic(:)=0.0 |
---|
| 345 | !----------------------------------------------------------------------- |
---|
| 346 | ! Choix du calendrier |
---|
| 347 | ! ------------------- |
---|
| 348 | |
---|
| 349 | ! calend = 'earth_365d' |
---|
| 350 | if (calend == 'earth_360d') then |
---|
| 351 | call ioconf_calendar('360d') |
---|
| 352 | write(*,*)'CALENDRIER CHOISI: Terrestre a 360 jours/an' |
---|
| 353 | else if (calend == 'earth_365d') then |
---|
| 354 | call ioconf_calendar('noleap') |
---|
| 355 | write(*,*)'CALENDRIER CHOISI: Terrestre a 365 jours/an' |
---|
| 356 | else if (calend == 'earth_366d') then |
---|
| 357 | call ioconf_calendar('all_leap') |
---|
| 358 | write(*,*)'CALENDRIER CHOISI: Terrestre bissextile' |
---|
| 359 | else if (calend == 'gregorian') then |
---|
| 360 | stop 'gregorian calend should not be used by normal user' |
---|
| 361 | call ioconf_calendar('gregorian') ! not to be used by normal users |
---|
| 362 | write(*,*)'CALENDRIER CHOISI: Gregorien' |
---|
| 363 | else |
---|
| 364 | write (*,*) 'ERROR : unknown calendar ', calend |
---|
| 365 | stop 'calend should be 360d,earth_365d,earth_366d,gregorian' |
---|
| 366 | endif |
---|
| 367 | !----------------------------------------------------------------------- |
---|
| 368 | ! |
---|
| 369 | !c Date : |
---|
| 370 | ! La date est supposee donnee sous la forme [annee, numero du jour dans |
---|
| 371 | ! l annee] ; l heure est donnee dans time_ini, lu dans lmdz1d.def. |
---|
| 372 | ! On appelle ymds2ju pour convertir [annee, jour] en [jour Julien]. |
---|
| 373 | ! Le numero du jour est dans "day". L heure est traitee separement. |
---|
| 374 | ! La date complete est dans "daytime" (l'unite est le jour). |
---|
| 375 | if (nday>0) then |
---|
| 376 | fnday=nday |
---|
| 377 | else |
---|
| 378 | fnday=-nday/float(day_step) |
---|
| 379 | endif |
---|
| 380 | print *,'fnday=',fnday |
---|
| 381 | ! start_time doit etre en FRACTION DE JOUR |
---|
| 382 | start_time=time_ini/24. |
---|
| 383 | |
---|
| 384 | annee_ref = anneeref |
---|
| 385 | mois = 1 |
---|
| 386 | day_ref = dayref |
---|
| 387 | heure = 0. |
---|
| 388 | itau_dyn = 0 |
---|
| 389 | itau_phy = 0 |
---|
| 390 | call ymds2ju(annee_ref,mois,day_ref,heure,day) |
---|
| 391 | day_ini = int(day) |
---|
| 392 | day_end = day_ini + int(fnday) |
---|
| 393 | |
---|
| 394 | ! Convert the initial date to Julian day |
---|
| 395 | day_ini_cas=day_deb |
---|
| 396 | print*,'time case',year_ini_cas,mth_ini_cas,day_ini_cas |
---|
| 397 | call ymds2ju & |
---|
| 398 | & (year_ini_cas,mth_ini_cas,day_ini_cas,heure_ini_cas*3600 & |
---|
| 399 | & ,day_ju_ini_cas) |
---|
| 400 | print*,'time case 2',day_ini_cas,day_ju_ini_cas |
---|
| 401 | daytime = day + heure_ini_cas/24. ! 1st day and initial time of the simulation |
---|
| 402 | |
---|
| 403 | ! Print out the actual date of the beginning of the simulation : |
---|
| 404 | call ju2ymds(daytime,year_print, month_print,day_print,sec_print) |
---|
| 405 | print *,' Time of beginning : ', & |
---|
| 406 | & year_print, month_print, day_print, sec_print |
---|
| 407 | |
---|
| 408 | !--------------------------------------------------------------------- |
---|
| 409 | ! Initialization of dimensions, geometry and initial state |
---|
| 410 | !--------------------------------------------------------------------- |
---|
| 411 | ! call init_phys_lmdz(1,1,llm,1,(/1/)) ! job now done via iniphysiq |
---|
| 412 | ! but we still need to initialize dimphy module (klon,klev,etc.) here. |
---|
| 413 | call init_dimphy1D(1,llm) |
---|
| 414 | call suphel |
---|
| 415 | call infotrac_init |
---|
| 416 | |
---|
| 417 | if (nqtot>nqmx) STOP 'Augmenter nqmx dans lmdz1d.F' |
---|
| 418 | allocate(q(llm,nqtot)) ; q(:,:)=0. |
---|
| 419 | allocate(dq(llm,nqtot)) |
---|
[3595] | 420 | allocate(d_q_vert_adv(llm,nqtot)) |
---|
[3541] | 421 | allocate(d_q_adv(llm,nqtot)) |
---|
| 422 | allocate(d_q_nudge(llm,nqtot)) |
---|
| 423 | ! allocate(d_th_adv(llm)) |
---|
| 424 | |
---|
| 425 | q(:,:) = 0. |
---|
| 426 | dq(:,:) = 0. |
---|
[3595] | 427 | d_q_vert_adv(:,:) = 0. |
---|
[3541] | 428 | d_q_adv(:,:) = 0. |
---|
| 429 | d_q_nudge(:,:) = 0. |
---|
| 430 | |
---|
| 431 | ! |
---|
| 432 | ! No ozone climatology need be read in this pre-initialization |
---|
| 433 | ! (phys_state_var_init is called again in physiq) |
---|
| 434 | read_climoz = 0 |
---|
| 435 | ! |
---|
| 436 | call phys_state_var_init(read_climoz) |
---|
| 437 | |
---|
| 438 | if (ngrid.ne.klon) then |
---|
| 439 | print*,'stop in inifis' |
---|
| 440 | print*,'Probleme de dimensions :' |
---|
| 441 | print*,'ngrid = ',ngrid |
---|
| 442 | print*,'klon = ',klon |
---|
| 443 | stop |
---|
| 444 | endif |
---|
| 445 | !!!===================================================================== |
---|
| 446 | !!! Feedback forcing values for Gateaux differentiation (al1) |
---|
| 447 | !!!===================================================================== |
---|
| 448 | !!! Surface Planck forcing bracketing call radiation |
---|
| 449 | !! surf_Planck = 0. |
---|
| 450 | !! surf_Conv = 0. |
---|
| 451 | !! write(*,*) 'Gateaux-dif Planck,Conv:',surf_Planck,surf_Conv |
---|
| 452 | !!! a mettre dans le lmdz1d.def ou autre |
---|
| 453 | !! |
---|
| 454 | !! |
---|
| 455 | qsol = qsolinp |
---|
| 456 | qsurf = fq_sat(tsurf,psurf/100.) |
---|
| 457 | day1= day_ini |
---|
| 458 | time=daytime-day |
---|
| 459 | ts_toga(1)=tsurf ! needed by read_tsurf1d.F |
---|
| 460 | rho(1)=psurf/(rd*tsurf*(1.+(rv/rd-1.)*qsurf)) |
---|
| 461 | |
---|
| 462 | ! |
---|
| 463 | !! mpl et jyg le 22/08/2012 : |
---|
| 464 | !! pour que les cas a flux de surface imposes marchent |
---|
| 465 | IF(.NOT.ok_flux_surf.or.max(abs(wtsurf),abs(wqsurf))>0.) THEN |
---|
| 466 | fsens=-wtsurf*rcpd*rho(1) |
---|
| 467 | flat=-wqsurf*rlvtt*rho(1) |
---|
| 468 | print *,'Flux: ok_flux wtsurf wqsurf',ok_flux_surf,wtsurf,wqsurf |
---|
| 469 | ENDIF |
---|
| 470 | print*,'Flux sol ',fsens,flat |
---|
| 471 | !! ok_flux_surf=.false. |
---|
| 472 | !! fsens=-wtsurf*rcpd*rho(1) |
---|
| 473 | !! flat=-wqsurf*rlvtt*rho(1) |
---|
| 474 | !!!! |
---|
| 475 | |
---|
| 476 | ! Vertical discretization and pressure levels at half and mid levels: |
---|
| 477 | |
---|
| 478 | pa = 5e4 |
---|
| 479 | !! preff= 1.01325e5 |
---|
| 480 | preff = psurf |
---|
| 481 | IF (ok_old_disvert) THEN |
---|
| 482 | call disvert0(pa,preff,ap,bp,dpres,presnivs,nivsigs,nivsig) |
---|
| 483 | print *,'On utilise disvert0' |
---|
| 484 | aps(1:llm)=0.5*(ap(1:llm)+ap(2:llm+1)) |
---|
| 485 | bps(1:llm)=0.5*(bp(1:llm)+bp(2:llm+1)) |
---|
| 486 | scaleheight=8. |
---|
| 487 | pseudoalt(1:llm)=-scaleheight*log(presnivs(1:llm)/preff) |
---|
| 488 | ELSE |
---|
| 489 | call disvert() |
---|
| 490 | print *,'On utilise disvert' |
---|
| 491 | ! Nouvelle version disvert permettant d imposer ap,bp (modif L.Guez) MPL 18092012 |
---|
| 492 | ! Dans ce cas, on lit ap,bp dans le fichier hybrid.txt |
---|
| 493 | ENDIF |
---|
| 494 | |
---|
| 495 | sig_s=presnivs/preff |
---|
| 496 | plev =ap+bp*psurf |
---|
| 497 | play = 0.5*(plev(1:llm)+plev(2:llm+1)) |
---|
| 498 | zlay=-rd*300.*log(play/psurf)/rg ! moved after reading profiles |
---|
| 499 | |
---|
| 500 | IF (forcing_type .eq. 59) THEN |
---|
| 501 | ! pour forcing_sandu, on cherche l'indice le plus proche de 700hpa#3000m |
---|
| 502 | write(*,*) '***********************' |
---|
| 503 | do l = 1, llm |
---|
| 504 | write(*,*) 'l,play(l),presnivs(l): ',l,play(l),presnivs(l) |
---|
| 505 | if (trouve_700 .and. play(l).le.70000) then |
---|
| 506 | llm700=l |
---|
| 507 | print *,'llm700,play=',llm700,play(l)/100. |
---|
| 508 | trouve_700= .false. |
---|
| 509 | endif |
---|
| 510 | enddo |
---|
| 511 | write(*,*) '***********************' |
---|
| 512 | ENDIF |
---|
| 513 | |
---|
| 514 | ! |
---|
| 515 | !===================================================================== |
---|
| 516 | ! EVENTUALLY, READ FORCING DATA : |
---|
| 517 | !===================================================================== |
---|
| 518 | |
---|
| 519 | #include "1D_read_forc_cases.h" |
---|
[3592] | 520 | print*,'A d_t_adv ',d_t_adv(1:20)*86400 |
---|
[3541] | 521 | |
---|
| 522 | if (forcing_GCM2SCM) then |
---|
| 523 | write (*,*) 'forcing_GCM2SCM not yet implemented' |
---|
| 524 | stop 'in initialization' |
---|
| 525 | endif ! forcing_GCM2SCM |
---|
| 526 | |
---|
| 527 | print*,'mxcalc=',mxcalc |
---|
| 528 | ! print*,'zlay=',zlay(mxcalc) |
---|
| 529 | print*,'play=',play(mxcalc) |
---|
| 530 | |
---|
| 531 | !Al1 pour SST forced, appell?? depuis ocean_forced_noice |
---|
| 532 | ts_cur = tsurf ! SST used in read_tsurf1d |
---|
| 533 | !===================================================================== |
---|
| 534 | ! Initialisation de la physique : |
---|
| 535 | !===================================================================== |
---|
| 536 | |
---|
| 537 | ! Rq: conf_phys.F90 lit tous les flags de physiq.def; conf_phys appele depuis physiq.F |
---|
| 538 | ! |
---|
| 539 | ! day_step, iphysiq lus dans gcm.def ci-dessus |
---|
| 540 | ! timestep: calcule ci-dessous from rday et day_step |
---|
| 541 | ! ngrid=1 |
---|
| 542 | ! llm: defini dans .../modipsl/modeles/LMDZ4/libf/grid/dimension |
---|
| 543 | ! rday: defini dans suphel.F (86400.) |
---|
| 544 | ! day_ini: lu dans run.def (dayref) |
---|
| 545 | ! rlat_rad,rlon-rad: transformes en radian de rlat,rlon lus dans lmdz1d.def (en degres) |
---|
| 546 | ! airefi,zcufi,zcvfi initialises au debut de ce programme |
---|
| 547 | ! rday,ra,rg,rd,rcpd declares dans YOMCST.h et calcules dans suphel.F |
---|
| 548 | day_step = float(nsplit_phys)*day_step/float(iphysiq) |
---|
| 549 | write (*,*) 'Time step divided by nsplit_phys (=',nsplit_phys,')' |
---|
| 550 | timestep =rday/day_step |
---|
| 551 | dtime_frcg = timestep |
---|
| 552 | ! |
---|
| 553 | zcufi=airefi |
---|
| 554 | zcvfi=airefi |
---|
| 555 | ! |
---|
| 556 | rlat_rad(1)=xlat*rpi/180. |
---|
| 557 | rlon_rad(1)=xlon*rpi/180. |
---|
| 558 | |
---|
| 559 | ! iniphysiq will call iniaqua who needs year_len from phys_cal_mod |
---|
| 560 | year_len_phys_cal_mod=year_len |
---|
| 561 | |
---|
| 562 | ! Ehouarn: iniphysiq requires arrays related to (3D) dynamics grid, |
---|
| 563 | ! e.g. for cell boundaries, which are meaningless in 1D; so pad these |
---|
| 564 | ! with '0.' when necessary |
---|
| 565 | call iniphysiq(iim,jjm,llm, & |
---|
| 566 | 1,comm_lmdz, & |
---|
| 567 | rday,day_ini,timestep, & |
---|
| 568 | (/rlat_rad(1),0./),(/0./), & |
---|
| 569 | (/0.,0./),(/rlon_rad(1),0./), & |
---|
| 570 | (/ (/airefi,0./),(/0.,0./) /), & |
---|
| 571 | (/zcufi,0.,0.,0./), & |
---|
| 572 | (/zcvfi,0./), & |
---|
| 573 | ra,rg,rd,rcpd,1) |
---|
| 574 | print*,'apres iniphysiq' |
---|
| 575 | |
---|
| 576 | ! 2 PARAMETRES QUI DEVRAIENT ETRE LUS DANS run.def MAIS NE LE SONT PAS ICI: |
---|
| 577 | co2_ppm= 330.0 |
---|
| 578 | solaire=1370.0 |
---|
| 579 | |
---|
| 580 | ! Ecriture du startphy avant le premier appel a la physique. |
---|
| 581 | ! On le met juste avant pour avoir acces a tous les champs |
---|
| 582 | |
---|
| 583 | if (ok_writedem) then |
---|
| 584 | |
---|
| 585 | !-------------------------------------------------------------------------- |
---|
| 586 | ! pbl_surface_init (called here) and pbl_surface_final (called by phyredem) |
---|
| 587 | ! need : qsol fder snow qsurf evap rugos agesno ftsoil |
---|
| 588 | !-------------------------------------------------------------------------- |
---|
| 589 | |
---|
| 590 | type_ocean = "force" |
---|
| 591 | run_off_lic_0(1) = restart_runoff |
---|
| 592 | call fonte_neige_init(run_off_lic_0) |
---|
| 593 | |
---|
| 594 | fder=0. |
---|
| 595 | snsrf(1,:)=snowmass ! masse de neige des sous surface |
---|
| 596 | qsurfsrf(1,:)=qsurf ! humidite de l'air des sous surface |
---|
| 597 | fevap=0. |
---|
| 598 | z0m(1,:)=rugos ! couverture de neige des sous surface |
---|
| 599 | z0h(1,:)=rugosh ! couverture de neige des sous surface |
---|
| 600 | agesno = xagesno |
---|
| 601 | tsoil(:,:,:)=tsurf |
---|
| 602 | !------ AMMA 2e run avec modele sol et rayonnement actif (MPL 23052012) |
---|
| 603 | ! tsoil(1,1,1)=299.18 |
---|
| 604 | ! tsoil(1,2,1)=300.08 |
---|
| 605 | ! tsoil(1,3,1)=301.88 |
---|
| 606 | ! tsoil(1,4,1)=305.48 |
---|
| 607 | ! tsoil(1,5,1)=308.00 |
---|
| 608 | ! tsoil(1,6,1)=308.00 |
---|
| 609 | ! tsoil(1,7,1)=308.00 |
---|
| 610 | ! tsoil(1,8,1)=308.00 |
---|
| 611 | ! tsoil(1,9,1)=308.00 |
---|
| 612 | ! tsoil(1,10,1)=308.00 |
---|
| 613 | ! tsoil(1,11,1)=308.00 |
---|
| 614 | !----------------------------------------------------------------------- |
---|
| 615 | call pbl_surface_init(fder, snsrf, qsurfsrf, tsoil) |
---|
| 616 | |
---|
| 617 | !------------------ prepare limit conditions for limit.nc ----------------- |
---|
| 618 | !-- Ocean force |
---|
| 619 | |
---|
| 620 | print*,'avant phyredem' |
---|
| 621 | pctsrf(1,:)=0. |
---|
| 622 | if (nat_surf.eq.0.) then |
---|
| 623 | pctsrf(1,is_oce)=1. |
---|
| 624 | pctsrf(1,is_ter)=0. |
---|
| 625 | pctsrf(1,is_lic)=0. |
---|
| 626 | pctsrf(1,is_sic)=0. |
---|
| 627 | else if (nat_surf .eq. 1) then |
---|
| 628 | pctsrf(1,is_oce)=0. |
---|
| 629 | pctsrf(1,is_ter)=1. |
---|
| 630 | pctsrf(1,is_lic)=0. |
---|
| 631 | pctsrf(1,is_sic)=0. |
---|
| 632 | else if (nat_surf .eq. 2) then |
---|
| 633 | pctsrf(1,is_oce)=0. |
---|
| 634 | pctsrf(1,is_ter)=0. |
---|
| 635 | pctsrf(1,is_lic)=1. |
---|
| 636 | pctsrf(1,is_sic)=0. |
---|
| 637 | else if (nat_surf .eq. 3) then |
---|
| 638 | pctsrf(1,is_oce)=0. |
---|
| 639 | pctsrf(1,is_ter)=0. |
---|
| 640 | pctsrf(1,is_lic)=0. |
---|
| 641 | pctsrf(1,is_sic)=1. |
---|
| 642 | |
---|
| 643 | end if |
---|
| 644 | |
---|
| 645 | |
---|
| 646 | print*,'nat_surf,pctsrf(1,is_oce),pctsrf(1,is_ter)',nat_surf & |
---|
| 647 | & ,pctsrf(1,is_oce),pctsrf(1,is_ter) |
---|
| 648 | |
---|
| 649 | zmasq=pctsrf(1,is_ter)+pctsrf(1,is_lic) |
---|
| 650 | zpic = zpicinp |
---|
| 651 | ftsol=tsurf |
---|
| 652 | nsw=6 ! on met le nb de bandes SW=6, pour initialiser |
---|
| 653 | ! 6 albedo, mais on peut quand meme tourner avec |
---|
| 654 | ! moins. Seules les 2 ou 4 premiers seront lus |
---|
| 655 | falb_dir=albedo |
---|
| 656 | falb_dif=albedo |
---|
| 657 | rugoro=rugos |
---|
| 658 | t_ancien(1,:)=temp(:) |
---|
| 659 | q_ancien(1,:)=q(:,1) |
---|
| 660 | ql_ancien = 0. |
---|
| 661 | qs_ancien = 0. |
---|
| 662 | prlw_ancien = 0. |
---|
| 663 | prsw_ancien = 0. |
---|
| 664 | prw_ancien = 0. |
---|
| 665 | !jyg< |
---|
| 666 | !! pbl_tke(:,:,:)=1.e-8 |
---|
| 667 | pbl_tke(:,:,:)=0. |
---|
| 668 | pbl_tke(:,2,:)=1.e-2 |
---|
| 669 | PRINT *, ' pbl_tke dans lmdz1d ' |
---|
| 670 | if (prt_level .ge. 5) then |
---|
| 671 | DO nsrf = 1,4 |
---|
| 672 | PRINT *,'pbl_tke(1,:,',nsrf,') ',pbl_tke(1,:,nsrf) |
---|
| 673 | ENDDO |
---|
| 674 | end if |
---|
| 675 | |
---|
| 676 | !>jyg |
---|
| 677 | |
---|
| 678 | rain_fall=0. |
---|
| 679 | snow_fall=0. |
---|
| 680 | solsw=0. |
---|
| 681 | sollw=0. |
---|
| 682 | sollwdown=rsigma*tsurf**4 |
---|
| 683 | radsol=0. |
---|
| 684 | rnebcon=0. |
---|
| 685 | ratqs=0. |
---|
| 686 | clwcon=0. |
---|
| 687 | zmax0 = 0. |
---|
| 688 | zmea=0. |
---|
| 689 | zstd=0. |
---|
| 690 | zsig=0. |
---|
| 691 | zgam=0. |
---|
| 692 | zval=0. |
---|
| 693 | zthe=0. |
---|
| 694 | sig1=0. |
---|
| 695 | w01=0. |
---|
| 696 | wake_cstar = 0. |
---|
| 697 | wake_deltaq = 0. |
---|
| 698 | wake_deltat = 0. |
---|
| 699 | wake_delta_pbl_TKE(:,:,:) = 0. |
---|
| 700 | delta_tsurf = 0. |
---|
| 701 | wake_fip = 0. |
---|
| 702 | wake_pe = 0. |
---|
| 703 | wake_s = 0. |
---|
| 704 | wake_dens = 0. |
---|
| 705 | ale_bl = 0. |
---|
| 706 | ale_bl_trig = 0. |
---|
| 707 | alp_bl = 0. |
---|
| 708 | IF (ALLOCATED(du_gwd_rando)) du_gwd_rando = 0. |
---|
| 709 | IF (ALLOCATED(du_gwd_front)) du_gwd_front = 0. |
---|
| 710 | entr_therm = 0. |
---|
| 711 | detr_therm = 0. |
---|
| 712 | f0 = 0. |
---|
| 713 | fm_therm = 0. |
---|
| 714 | u_ancien(1,:)=u(:) |
---|
| 715 | v_ancien(1,:)=v(:) |
---|
| 716 | |
---|
[3592] | 717 | u10m=0. |
---|
| 718 | v10m=0. |
---|
| 719 | ale_wake=0. |
---|
| 720 | ale_bl_stat=0. |
---|
| 721 | |
---|
[3541] | 722 | !------------------------------------------------------------------------ |
---|
| 723 | ! Make file containing restart for the physics (startphy.nc) |
---|
| 724 | ! |
---|
| 725 | ! NB: List of the variables to be written by phyredem (via put_field): |
---|
| 726 | ! rlon,rlat,zmasq,pctsrf(:,is_ter),pctsrf(:,is_lic),pctsrf(:,is_oce) |
---|
| 727 | ! pctsrf(:,is_sic),ftsol(:,nsrf),tsoil(:,isoil,nsrf),qsurf(:,nsrf) |
---|
| 728 | ! qsol,falb_dir(:,nsrf),falb_dif(:,nsrf),evap(:,nsrf),snow(:,nsrf) |
---|
| 729 | ! radsol,solsw,sollw, sollwdown,fder,rain_fall,snow_fall,frugs(:,nsrf) |
---|
| 730 | ! agesno(:,nsrf),zmea,zstd,zsig,zgam,zthe,zpic,zval,rugoro |
---|
| 731 | ! t_ancien,q_ancien,,frugs(:,is_oce),clwcon(:,1),rnebcon(:,1),ratqs(:,1) |
---|
| 732 | ! run_off_lic_0,pbl_tke(:,1:klev,nsrf), zmax0,f0,sig1,w01 |
---|
| 733 | ! wake_deltat,wake_deltaq,wake_s,wake_dens,wake_cstar, |
---|
| 734 | ! wake_fip,wake_delta_pbl_tke(:,1:klev,nsrf) |
---|
| 735 | ! |
---|
| 736 | ! NB2: The content of the startphy.nc file depends on some flags defined in |
---|
| 737 | ! the ".def" files. However, since conf_phys is not called in lmdz1d.F90, these flags have |
---|
| 738 | ! to be set at some arbitratry convenient values. |
---|
| 739 | !------------------------------------------------------------------------ |
---|
| 740 | !Al1 =============== restart option ========================== |
---|
| 741 | if (.not.restart) then |
---|
| 742 | iflag_pbl = 5 |
---|
| 743 | call phyredem ("startphy.nc") |
---|
| 744 | else |
---|
| 745 | ! (desallocations) |
---|
| 746 | print*,'callin surf final' |
---|
| 747 | call pbl_surface_final( fder, snsrf, qsurfsrf, tsoil) |
---|
| 748 | print*,'after surf final' |
---|
| 749 | CALL fonte_neige_final(run_off_lic_0) |
---|
| 750 | endif |
---|
| 751 | |
---|
| 752 | ok_writedem=.false. |
---|
| 753 | print*,'apres phyredem' |
---|
| 754 | |
---|
| 755 | endif ! ok_writedem |
---|
| 756 | |
---|
| 757 | !------------------------------------------------------------------------ |
---|
| 758 | ! Make file containing boundary conditions (limit.nc) **Al1->restartdyn*** |
---|
| 759 | ! -------------------------------------------------- |
---|
| 760 | ! NB: List of the variables to be written in limit.nc |
---|
| 761 | ! (by writelim.F, subroutine of 1DUTILS.h): |
---|
| 762 | ! phy_nat,phy_alb,phy_sst,phy_bil,phy_rug,phy_ice, |
---|
| 763 | ! phy_fter,phy_foce,phy_flic,phy_fsic) |
---|
| 764 | !------------------------------------------------------------------------ |
---|
| 765 | do i=1,year_len |
---|
| 766 | phy_nat(i) = nat_surf |
---|
| 767 | phy_alb(i) = albedo |
---|
| 768 | phy_sst(i) = tsurf ! read_tsurf1d will be used instead |
---|
| 769 | phy_rug(i) = rugos |
---|
| 770 | phy_fter(i) = pctsrf(1,is_ter) |
---|
| 771 | phy_foce(i) = pctsrf(1,is_oce) |
---|
| 772 | phy_fsic(i) = pctsrf(1,is_sic) |
---|
| 773 | phy_flic(i) = pctsrf(1,is_lic) |
---|
| 774 | enddo |
---|
| 775 | |
---|
| 776 | ! fabrication de limit.nc |
---|
| 777 | call writelim (1,phy_nat,phy_alb,phy_sst,phy_bil,phy_rug, & |
---|
| 778 | & phy_ice,phy_fter,phy_foce,phy_flic,phy_fsic) |
---|
| 779 | |
---|
| 780 | |
---|
| 781 | call phys_state_var_end |
---|
| 782 | !Al1 |
---|
| 783 | if (restart) then |
---|
| 784 | print*,'call to restart dyn 1d' |
---|
| 785 | Call dyn1deta0("start1dyn.nc",plev,play,phi,phis,presnivs, & |
---|
| 786 | & u,v,temp,q,omega2) |
---|
| 787 | |
---|
| 788 | print*,'fnday,annee_ref,day_ref,day_ini', & |
---|
| 789 | & fnday,annee_ref,day_ref,day_ini |
---|
| 790 | !** call ymds2ju(annee_ref,mois,day_ini,heure,day) |
---|
| 791 | day = day_ini |
---|
| 792 | day_end = day_ini + nday |
---|
| 793 | daytime = day + time_ini/24. ! 1st day and initial time of the simulation |
---|
| 794 | |
---|
| 795 | ! Print out the actual date of the beginning of the simulation : |
---|
| 796 | call ju2ymds(daytime, an, mois, jour, heure) |
---|
| 797 | print *,' Time of beginning : y m d h',an, mois,jour,heure/3600. |
---|
| 798 | |
---|
| 799 | day = int(daytime) |
---|
| 800 | time=daytime-day |
---|
| 801 | |
---|
| 802 | print*,'****** intialised fields from restart1dyn *******' |
---|
| 803 | print*,'plev,play,phi,phis,presnivs,u,v,temp,q,omega2' |
---|
| 804 | print*,'temp(1),q(1,1),u(1),v(1),plev(1),phis :' |
---|
| 805 | print*,temp(1),q(1,1),u(1),v(1),plev(1),phis |
---|
| 806 | ! raz for safety |
---|
| 807 | do l=1,llm |
---|
[3595] | 808 | d_q_vert_adv(l,1) = 0. |
---|
[3541] | 809 | enddo |
---|
| 810 | endif |
---|
| 811 | !Al1 ================ end restart ================================= |
---|
| 812 | IF (ecrit_slab_oc.eq.1) then |
---|
| 813 | open(97,file='div_slab.dat',STATUS='UNKNOWN') |
---|
| 814 | elseif (ecrit_slab_oc.eq.0) then |
---|
| 815 | open(97,file='div_slab.dat',STATUS='OLD') |
---|
| 816 | endif |
---|
| 817 | ! |
---|
| 818 | !===================================================================== |
---|
[3595] | 819 | #ifdef OUTPUT_PHYS_SCM |
---|
[3541] | 820 | CALL iophys_ini |
---|
[3595] | 821 | #endif |
---|
[3541] | 822 | ! START OF THE TEMPORAL LOOP : |
---|
| 823 | !===================================================================== |
---|
| 824 | |
---|
| 825 | it_end = nint(fnday*day_step) |
---|
| 826 | !test JLD it_end = 10 |
---|
| 827 | do while(it.le.it_end) |
---|
| 828 | |
---|
| 829 | if (prt_level.ge.1) then |
---|
| 830 | print*,'XXXXXXXXXXXXXXXXXXX ITAP,day,time=', & |
---|
| 831 | & it,day,time,it_end,day_step |
---|
| 832 | print*,'PAS DE TEMPS ',timestep |
---|
| 833 | endif |
---|
| 834 | !Al1 demande de restartphy.nc |
---|
| 835 | if (it.eq.it_end) lastcall=.True. |
---|
| 836 | |
---|
| 837 | !--------------------------------------------------------------------- |
---|
| 838 | ! Interpolation of forcings in time and onto model levels |
---|
| 839 | !--------------------------------------------------------------------- |
---|
| 840 | |
---|
| 841 | #include "1D_interp_cases.h" |
---|
| 842 | |
---|
| 843 | !--------------------------------------------------------------------- |
---|
| 844 | ! Geopotential : |
---|
| 845 | !--------------------------------------------------------------------- |
---|
| 846 | |
---|
| 847 | phi(1)=RD*temp(1)*(plev(1)-play(1))/(.5*(plev(1)+play(1))) |
---|
| 848 | do l = 1, llm-1 |
---|
| 849 | phi(l+1)=phi(l)+RD*(temp(l)+temp(l+1))* & |
---|
| 850 | & (play(l)-play(l+1))/(play(l)+play(l+1)) |
---|
| 851 | enddo |
---|
| 852 | |
---|
| 853 | !--------------------------------------------------------------------- |
---|
[3595] | 854 | ! Vertical advection |
---|
| 855 | !--------------------------------------------------------------------- |
---|
| 856 | |
---|
| 857 | IF ( forc_w+forc_omega > 0 ) THEN |
---|
| 858 | |
---|
| 859 | IF ( forc_w == 1 ) THEN |
---|
| 860 | w_adv=w_mod_cas |
---|
| 861 | z_adv=phi/RG |
---|
| 862 | ELSE |
---|
| 863 | w_adv=omega |
---|
| 864 | z_adv=play |
---|
| 865 | ENDIF |
---|
| 866 | |
---|
| 867 | teta=temp*(pzero/play)**rkappa |
---|
| 868 | do l=2,llm-1 |
---|
[3693] | 869 | ! vertical tendencies computed as d X / d t = -W d X / d z |
---|
[3595] | 870 | d_u_vert_adv(l)=-w_adv(l)*(u(l+1)-u(l-1))/(z_adv(l+1)-z_adv(l-1)) |
---|
| 871 | d_v_vert_adv(l)=-w_adv(l)*(v(l+1)-v(l-1))/(z_adv(l+1)-z_adv(l-1)) |
---|
[3693] | 872 | ! d theta / dt = -W d theta / d z, transformed into d temp / d t dividing by (pzero/play(l))**rkappa |
---|
| 873 | 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 |
---|
[3595] | 874 | 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)) |
---|
| 875 | enddo |
---|
[3693] | 876 | d_u_adv(:)=d_u_adv(:)+d_u_vert_adv(:) |
---|
| 877 | d_v_adv(:)=d_v_adv(:)+d_v_vert_adv(:) |
---|
[3595] | 878 | d_t_adv(:)=d_t_adv(:)+d_t_vert_adv(:) |
---|
| 879 | d_q_adv(:,1)=d_q_adv(:,1)+d_q_vert_adv(:,1) |
---|
| 880 | |
---|
| 881 | print*,'OMEGA ',w_adv(10),z_adv(10) |
---|
| 882 | |
---|
| 883 | ENDIF |
---|
| 884 | |
---|
| 885 | !--------------------------------------------------------------------- |
---|
[3541] | 886 | ! Listing output for debug prt_level>=1 |
---|
| 887 | !--------------------------------------------------------------------- |
---|
| 888 | if (prt_level>=1) then |
---|
| 889 | print *,' avant physiq : -------- day time ',day,time |
---|
| 890 | write(*,*) 'firstcall,lastcall,phis', & |
---|
| 891 | & firstcall,lastcall,phis |
---|
| 892 | end if |
---|
| 893 | if (prt_level>=5) then |
---|
| 894 | write(*,'(a10,2a4,4a13)') 'BEFOR1 IT=','it','l', & |
---|
| 895 | & 'presniv','plev','play','phi' |
---|
| 896 | write(*,'(a10,2i4,4f13.2)') ('BEFOR1 IT= ',it,l, & |
---|
| 897 | & presnivs(l),plev(l),play(l),phi(l),l=1,llm) |
---|
| 898 | write(*,'(a11,2a4,a11,6a8)') 'BEFOR2','it','l', & |
---|
| 899 | & 'presniv','u','v','temp','q1','q2','omega2' |
---|
| 900 | write(*,'(a11,2i4,f11.2,5f8.2,e10.2)') ('BEFOR2 IT= ',it,l, & |
---|
| 901 | & presnivs(l),u(l),v(l),temp(l),q(l,1),q(l,2),omega2(l),l=1,llm) |
---|
| 902 | endif |
---|
| 903 | |
---|
| 904 | !--------------------------------------------------------------------- |
---|
| 905 | ! Call physiq : |
---|
| 906 | !--------------------------------------------------------------------- |
---|
| 907 | call physiq(ngrid,llm, & |
---|
| 908 | firstcall,lastcall,timestep, & |
---|
| 909 | plev,play,phi,phis,presnivs, & |
---|
| 910 | u,v, rot, temp,q,omega2, & |
---|
| 911 | du_phys,dv_phys,dt_phys,dq,dpsrf) |
---|
| 912 | firstcall=.false. |
---|
| 913 | |
---|
| 914 | !--------------------------------------------------------------------- |
---|
| 915 | ! Listing output for debug |
---|
| 916 | !--------------------------------------------------------------------- |
---|
| 917 | if (prt_level>=5) then |
---|
| 918 | write(*,'(a11,2a4,4a13)') 'AFTER1 IT=','it','l', & |
---|
| 919 | & 'presniv','plev','play','phi' |
---|
| 920 | write(*,'(a11,2i4,4f13.2)') ('AFTER1 it= ',it,l, & |
---|
| 921 | & presnivs(l),plev(l),play(l),phi(l),l=1,llm) |
---|
| 922 | write(*,'(a11,2a4,a11,6a8)') 'AFTER2','it','l', & |
---|
| 923 | & 'presniv','u','v','temp','q1','q2','omega2' |
---|
| 924 | write(*,'(a11,2i4,f11.2,5f8.2,e10.2)') ('AFTER2 it= ',it,l, & |
---|
| 925 | & presnivs(l),u(l),v(l),temp(l),q(l,1),q(l,2),omega2(l),l=1,llm) |
---|
| 926 | write(*,'(a11,2a4,a11,5a8)') 'AFTER3','it','l', & |
---|
| 927 | & 'presniv','du_phys','dv_phys','dt_phys','dq1','dq2' |
---|
| 928 | write(*,'(a11,2i4,f11.2,5f8.2)') ('AFTER3 it= ',it,l, & |
---|
| 929 | & presnivs(l),86400*du_phys(l),86400*dv_phys(l), & |
---|
| 930 | & 86400*dt_phys(l),86400*dq(l,1),dq(l,2),l=1,llm) |
---|
| 931 | write(*,*) 'dpsrf',dpsrf |
---|
| 932 | endif |
---|
| 933 | !--------------------------------------------------------------------- |
---|
| 934 | ! Add physical tendencies : |
---|
| 935 | !--------------------------------------------------------------------- |
---|
| 936 | |
---|
| 937 | fcoriolis=2.*sin(rpi*xlat/180.)*romega |
---|
| 938 | |
---|
| 939 | IF (prt_level >= 5) print*, 'fcoriolis, xlat,mxcalc ', & |
---|
| 940 | fcoriolis, xlat,mxcalc |
---|
| 941 | |
---|
[3693] | 942 | !--------------------------------------------------------------------- |
---|
| 943 | ! Geostrophic forcing |
---|
| 944 | !--------------------------------------------------------------------- |
---|
[3541] | 945 | |
---|
[3693] | 946 | IF ( forc_geo == 0 ) THEN |
---|
| 947 | d_u_age(1:mxcalc)=0. |
---|
| 948 | d_v_age(1:mxcalc)=0. |
---|
| 949 | ELSE |
---|
[3541] | 950 | sfdt = sin(0.5*fcoriolis*timestep) |
---|
| 951 | cfdt = cos(0.5*fcoriolis*timestep) |
---|
| 952 | ! print *,'fcoriolis,sfdt,cfdt,timestep',fcoriolis,sfdt,cfdt,timestep |
---|
| 953 | ! |
---|
[3693] | 954 | d_u_age(1:mxcalc)= -2.*sfdt/timestep* & |
---|
[3541] | 955 | & (sfdt*(u(1:mxcalc)-ug(1:mxcalc)) - & |
---|
| 956 | & cfdt*(v(1:mxcalc)-vg(1:mxcalc)) ) |
---|
| 957 | !! : fcoriolis*(v(1:mxcalc)-vg(1:mxcalc)) |
---|
| 958 | ! |
---|
[3693] | 959 | d_v_age(1:mxcalc)= -2.*sfdt/timestep* & |
---|
[3541] | 960 | & (cfdt*(u(1:mxcalc)-ug(1:mxcalc)) + & |
---|
| 961 | & sfdt*(v(1:mxcalc)-vg(1:mxcalc)) ) |
---|
| 962 | !! : -fcoriolis*(u(1:mxcalc)-ug(1:mxcalc)) |
---|
[3693] | 963 | ENDIF |
---|
[3541] | 964 | ! |
---|
[3693] | 965 | !--------------------------------------------------------------------- |
---|
[3541] | 966 | ! Nudging |
---|
[3693] | 967 | !--------------------------------------------------------------------- |
---|
[3541] | 968 | d_t_nudge(:) = 0. |
---|
| 969 | d_u_nudge(:) = 0. |
---|
| 970 | d_v_nudge(:) = 0. |
---|
[3593] | 971 | d_q_nudge(:,:) = 0. |
---|
[3594] | 972 | DO l=1,llm |
---|
| 973 | IF ( play(l) < p_nudging_u .AND. nint(nudging_u) /= 0 ) & |
---|
| 974 | & d_u_nudge(l)=(u_nudg_mod_cas(l)-u(l))/nudging_u |
---|
| 975 | IF ( play(l) < p_nudging_v .AND. nint(nudging_v) /= 0 ) & |
---|
| 976 | & d_v_nudge(l)=(v_nudg_mod_cas(l)-v(l))/nudging_v |
---|
[3686] | 977 | IF ( play(l) < p_nudging_t .AND. nint(nudging_t) /= 0 ) & |
---|
| 978 | & d_t_nudge(l)=(temp_nudg_mod_cas(l)-temp(l))/nudging_t |
---|
[3594] | 979 | IF ( play(l) < p_nudging_qv .AND. nint(nudging_qv) /= 0 ) & |
---|
| 980 | & d_q_nudge(l,1)=(qv_nudg_mod_cas(l)-q(l,1))/nudging_qv |
---|
| 981 | ENDDO |
---|
[3595] | 982 | |
---|
[3693] | 983 | !--------------------------------------------------------------------- |
---|
| 984 | ! Optional outputs |
---|
| 985 | !--------------------------------------------------------------------- |
---|
[3595] | 986 | #ifdef OUTPUT_PHYS_SCM |
---|
| 987 | CALL iophys_ecrit('w_adv',klev,'w_adv','K/day',w_adv) |
---|
| 988 | CALL iophys_ecrit('z_adv',klev,'z_adv','K/day',z_adv) |
---|
| 989 | CALL iophys_ecrit('dtadv',klev,'dtadv','K/day',86400*d_t_adv) |
---|
| 990 | CALL iophys_ecrit('dtdyn',klev,'dtdyn','K/day',86400*d_t_vert_adv) |
---|
[3594] | 991 | CALL iophys_ecrit('qv',klev,'qv','g/kg',1000*q(:,1)) |
---|
| 992 | CALL iophys_ecrit('qvnud',klev,'qvnud','g/kg',1000*u_nudg_mod_cas) |
---|
| 993 | CALL iophys_ecrit('u',klev,'u','m/s',u) |
---|
| 994 | CALL iophys_ecrit('unud',klev,'unud','m/s',u_nudg_mod_cas) |
---|
| 995 | CALL iophys_ecrit('v',klev,'v','m/s',v) |
---|
| 996 | CALL iophys_ecrit('vnud',klev,'vnud','m/s',v_nudg_mod_cas) |
---|
[3595] | 997 | CALL iophys_ecrit('temp',klev,'temp','K',temp) |
---|
| 998 | CALL iophys_ecrit('tempnud',klev,'temp_nudg_mod_cas','K',temp_nudg_mod_cas) |
---|
[3594] | 999 | CALL iophys_ecrit('dtnud',klev,'dtnud','K/day',86400*d_t_nudge) |
---|
| 1000 | CALL iophys_ecrit('dqnud',klev,'dqnud','K/day',1000*86400*d_q_nudge(:,1)) |
---|
[3595] | 1001 | #endif |
---|
[3593] | 1002 | |
---|
[3541] | 1003 | IF (flag_inhib_forcing == 0) then ! if tendency of forcings should be added |
---|
| 1004 | |
---|
| 1005 | u(1:mxcalc)=u(1:mxcalc) + timestep*( & |
---|
| 1006 | & du_phys(1:mxcalc) & |
---|
[3693] | 1007 | & +d_u_age(1:mxcalc)+d_u_adv(1:mxcalc) & |
---|
[3541] | 1008 | & +d_u_nudge(1:mxcalc) ) |
---|
| 1009 | v(1:mxcalc)=v(1:mxcalc) + timestep*( & |
---|
| 1010 | & dv_phys(1:mxcalc) & |
---|
[3693] | 1011 | & +d_v_age(1:mxcalc)+d_v_adv(1:mxcalc) & |
---|
[3541] | 1012 | & +d_v_nudge(1:mxcalc) ) |
---|
| 1013 | q(1:mxcalc,:)=q(1:mxcalc,:)+timestep*( & |
---|
| 1014 | & dq(1:mxcalc,:) & |
---|
| 1015 | & +d_q_adv(1:mxcalc,:) & |
---|
| 1016 | & +d_q_nudge(1:mxcalc,:) ) |
---|
| 1017 | |
---|
| 1018 | if (prt_level.ge.3) then |
---|
| 1019 | print *, & |
---|
| 1020 | & 'physiq-> temp(1),dt_phys(1),d_t_adv(1),dt_cooling(1) ', & |
---|
| 1021 | & temp(1),dt_phys(1),d_t_adv(1),dt_cooling(1) |
---|
| 1022 | print* ,'dv_phys=',dv_phys |
---|
[3693] | 1023 | print* ,'d_v_age=',d_v_age |
---|
| 1024 | print* ,'d_v_adv=',d_v_adv |
---|
[3541] | 1025 | print* ,'d_v_nudge=',d_v_nudge |
---|
| 1026 | print*, v |
---|
| 1027 | print*, vg |
---|
| 1028 | endif |
---|
| 1029 | |
---|
| 1030 | temp(1:mxcalc)=temp(1:mxcalc)+timestep*( & |
---|
| 1031 | & dt_phys(1:mxcalc) & |
---|
| 1032 | & +d_t_adv(1:mxcalc) & |
---|
| 1033 | & +d_t_nudge(1:mxcalc) & |
---|
| 1034 | & +dt_cooling(1:mxcalc)) ! Taux de chauffage ou refroid. |
---|
| 1035 | |
---|
| 1036 | |
---|
[3592] | 1037 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 1038 | !! CONSERVE EN ATTENDANT QUE LE CAS EN QUESTION FONCTIONNE EN STD !! |
---|
| 1039 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 1040 | ! endif ! forcing_sandu or forcing_astex |
---|
| 1041 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 1042 | |
---|
[3541] | 1043 | teta=temp*(pzero/play)**rkappa |
---|
| 1044 | ! |
---|
| 1045 | !--------------------------------------------------------------------- |
---|
| 1046 | ! Nudge soil temperature if requested |
---|
| 1047 | !--------------------------------------------------------------------- |
---|
| 1048 | |
---|
| 1049 | IF (nudge_tsoil .AND. .NOT. lastcall) THEN |
---|
| 1050 | ftsoil(1,isoil_nudge,:) = ftsoil(1,isoil_nudge,:) & |
---|
| 1051 | & -timestep/tau_soil_nudge*(ftsoil(1,isoil_nudge,:)-Tsoil_nudge) |
---|
| 1052 | ENDIF |
---|
| 1053 | |
---|
| 1054 | !--------------------------------------------------------------------- |
---|
| 1055 | ! Add large-scale tendencies (advection, etc) : |
---|
| 1056 | !--------------------------------------------------------------------- |
---|
| 1057 | |
---|
| 1058 | !cc nrlmd |
---|
| 1059 | !cc tmpvar=teta |
---|
| 1060 | !cc call advect_vert(llm,omega,timestep,tmpvar,plev) |
---|
| 1061 | !cc |
---|
| 1062 | !cc teta(1:mxcalc)=tmpvar(1:mxcalc) |
---|
| 1063 | !cc tmpvar(:)=q(:,1) |
---|
| 1064 | !cc call advect_vert(llm,omega,timestep,tmpvar,plev) |
---|
| 1065 | !cc q(1:mxcalc,1)=tmpvar(1:mxcalc) |
---|
| 1066 | !cc tmpvar(:)=q(:,2) |
---|
| 1067 | !cc call advect_vert(llm,omega,timestep,tmpvar,plev) |
---|
| 1068 | !cc q(1:mxcalc,2)=tmpvar(1:mxcalc) |
---|
| 1069 | |
---|
| 1070 | END IF ! end if tendency of tendency should be added |
---|
| 1071 | |
---|
| 1072 | !--------------------------------------------------------------------- |
---|
| 1073 | ! Air temperature : |
---|
| 1074 | !--------------------------------------------------------------------- |
---|
| 1075 | if (lastcall) then |
---|
| 1076 | print*,'Pas de temps final ',it |
---|
| 1077 | call ju2ymds(daytime, an, mois, jour, heure) |
---|
| 1078 | print*,'a la date : a m j h',an, mois, jour ,heure/3600. |
---|
| 1079 | endif |
---|
| 1080 | |
---|
| 1081 | ! incremente day time |
---|
| 1082 | ! print*,'daytime bef',daytime,1./day_step |
---|
| 1083 | daytime = daytime+1./day_step |
---|
| 1084 | !Al1dbg |
---|
| 1085 | day = int(daytime+0.1/day_step) |
---|
| 1086 | ! time = max(daytime-day,0.0) |
---|
| 1087 | !Al1&jyg: correction de bug |
---|
| 1088 | !cc time = real(mod(it,day_step))/day_step |
---|
| 1089 | time = time_ini/24.+real(mod(it,day_step))/day_step |
---|
| 1090 | ! print*,'daytime nxt time',daytime,time |
---|
| 1091 | it=it+1 |
---|
| 1092 | |
---|
| 1093 | enddo |
---|
| 1094 | |
---|
| 1095 | !Al1 |
---|
| 1096 | if (ecrit_slab_oc.ne.-1) close(97) |
---|
| 1097 | |
---|
| 1098 | !Al1 Call to 1D equivalent of dynredem (an,mois,jour,heure ?) |
---|
| 1099 | ! ------------------------------------- |
---|
| 1100 | call dyn1dredem("restart1dyn.nc", & |
---|
| 1101 | & plev,play,phi,phis,presnivs, & |
---|
| 1102 | & u,v,temp,q,omega2) |
---|
| 1103 | |
---|
| 1104 | CALL abort_gcm ('lmdz1d ','The End ',0) |
---|
| 1105 | |
---|
| 1106 | END SUBROUTINE scm |
---|