Changeset 3035 for trunk/LMDZ.VENUS/libf
- Timestamp:
- Sep 7, 2023, 3:47:46 PM (16 months ago)
- Location:
- trunk/LMDZ.VENUS/libf/phyvenus
- Files:
-
- 2 added
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.VENUS/libf/phyvenus/clesphys.h
r2836 r3035 11 11 LOGICAL ok_orodr,ok_orolf,ok_gw_nonoro 12 12 LOGICAL ok_kzmin,tuneupperatm, ok_ionchem, ok_jonline 13 LOGICAL ok_iondiff 13 14 LOGICAL callnlte,callnirco2,callthermos 14 15 LOGICAL ok_cloud, ok_chem, reinit_trac, ok_sedim 15 16 LOGICAL ok_clmain, physideal, startphy_file 16 INTEGER nbapp_rad, nbapp_ch im, iflag_con, iflag_ajs17 INTEGER nbapp_rad, nbapp_chem, iflag_con, iflag_ajs 17 18 INTEGER lev_histins, lev_histday, lev_histmth 18 19 INTEGER tr_scheme, cl_scheme … … 29 30 & tuneupperatm,callnlte,callnirco2,callthermos, & 30 31 & ok_cloud, ok_chem, reinit_trac, ok_sedim, & 31 & ok_clmain, physideal, startphy_file, ok_ionchem, ok_jonline 32 & ok_clmain, physideal, startphy_file, ok_ionchem, ok_jonline, & 33 & ok_iondiff 32 34 33 COMMON/clesphys_i/ nbapp_rad, nbapp_ch im, &35 COMMON/clesphys_i/ nbapp_rad, nbapp_chem, & 34 36 & iflag_con, iflag_ajs, & 35 37 & lev_histins, lev_histday, lev_histmth, tr_scheme, & -
trunk/LMDZ.VENUS/libf/phyvenus/conf_phys.F90
r2969 r3035 92 92 call getin('nbapp_rad',nbapp_rad) 93 93 print*,"nbapp_rad",nbapp_rad 94 !Config Key = nbapp_ch im94 !Config Key = nbapp_chem 95 95 !Config Desc = Frequence d'appel a la chimie 96 !Config Def = 196 !Config Def = 24000 97 97 !Config Help = Nombre d'appels des routines de chimie 98 98 !Config par jour. 99 nbapp_ch im = 1100 call getin('nbapp_ch im',nbapp_chim)99 nbapp_chem = 24000 100 call getin('nbapp_chem',nbapp_chem) 101 101 102 102 !Config Key = iflag_con … … 481 481 !Config Desc = 482 482 !Config Def = .false. 483 !Config Help = 483 !Config Help = Activate thermosphere tuning for oxygen compensation 484 484 ! 485 485 tuneupperatm = .false. … … 499 499 !Config Desc = 500 500 !Config Def = .false. 501 !Config Help = 501 !Config Help = Activate the ionosphere chemistry 502 502 ! 503 503 ok_ionchem = .false. … … 514 514 515 515 ! 516 !Config Key = ok_iondiff 517 !Config Desc = 518 !Config Def = .false. 519 !Config Help = Activate the ambipolar ion diffusion in physic 520 ! 521 ok_iondiff = .false. 522 call getin('ok_iondiff',ok_iondiff) 523 524 if ((ok_ionchem.eqv..false.).and.(.true..eqv.ok_iondiff)) then 525 write(*,*) "Error incoherent flags :" 526 write(*,*) "ok_ionchem=",ok_ionchem," / ok_iondiff=",ok_iondiff 527 write(*,*) "If you include ions diffusion," 528 write(*,*) " you need also ok_ionchem==.true." 529 write(*,*) "If you do not include ions diffusion," 530 write(*,*) " you need ok_iondiff==.false." 531 write(*,*) "Check physiq.def" 532 stop 533 endif 534 535 ! 516 536 ! 517 537 !Config Key = … … 535 555 write(lunout,*)' ok_gw_nonoro = ', ok_gw_nonoro 536 556 write(lunout,*)' nbapp_rad = ', nbapp_rad 537 write(lunout,*)' nbapp_ch im = ', nbapp_chim557 write(lunout,*)' nbapp_chem = ', nbapp_chem 538 558 write(lunout,*)' iflag_con = ', iflag_con 539 559 write(lunout,*)' Sortie journaliere = ', ok_journe … … 576 596 write(lunout,*)' ok_jonline = ',ok_jonline 577 597 write(lunout,*)' ok_ionchem = ',ok_ionchem 578 598 write(lunout,*)' ok_iondiff = ',ok_iondiff 599 579 600 end subroutine conf_phys 580 601 -
trunk/LMDZ.VENUS/libf/phyvenus/nonoro_gwd_ran_mod.F90
r2836 r3035 264 264 zp(jw, ii) = (sign(1., 0.5 - ran_num_1) & 265 265 + 1.) * RPI / 2. 266 ! zp(jw, ii) = atan2(4.*vh(ii, launch),uh(ii, launch)) & 267 ! + (sign(1., 0.5 - ran_num_1) + 1.) * RPI / 2. 266 268 ! --------- TRY 00 ----------------------- 267 269 !IF((abs(latitude_deg(ii)) .le. 55.)) THEN … … 301 303 sqrt(3.) / sqrt(na * 1.) 302 304 end DO 305 !cpha = cpha - uh(ii, launch) 303 306 IF (cpha < 0.) THEN 304 307 cpha = - 1. * cpha 305 308 zp (jw, ii) = zp(jw, ii) + RPI 306 309 ENDIF 310 !IF (cpha < 1.) THEN 311 ! cpha = 1. 312 !ENDIF 307 313 ! otherwise, with the computation below, we get a uniform distribution between cmin and cmax. 308 314 ! ran_num_3 = mod(tt(ii, jw)**2, 1.) -
trunk/LMDZ.VENUS/libf/phyvenus/photochemistry_venus.F90
r2949 r3035 2507 2507 !-- TuneE 2508 2508 ! VCD 2.0 tuning 2509 v_phot(65: 90,4) = v_phot(65:90,4)*10. ! CO2 + hv ==> O(1D) + CO2509 v_phot(65:nz,4) = v_phot(65:nz,4)*10. ! CO2 + hv ==> O(1D) + CO 2510 2510 !-- 2511 2511 ! v_phot(:,4) = v_phot(:,4)*10. … … 4590 4590 v_phot(65:nz,4) = v_phot(65:nz,4)*8. ! CO2 + hv ==> O(1D) + CO 4591 4591 !-- 4592 ! TEST TUNEH 4593 ! v_phot(62:74,4) = v_phot(62:74,4)*6.5 ! CO2 + hv ==> O(1D) + CO 4594 !-- 4592 4595 ! v_phot(:,4) = v_phot(:,4)*10. 4593 4596 !do ij=3,4 -
trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F
r2836 r3035 298 298 299 299 integer :: chempas 300 real :: nbapp_chem,zctime300 real :: zctime 301 301 302 302 ! sedimentation … … 340 340 c Tendencies due to molecular diffusion 341 341 real d_q_moldif(klon,klev,nqmax) 342 343 c Tendencies due to ambipolar ion diffusion 344 real d_q_iondif(klon,klev,nqmax) 342 345 343 346 c … … 671 674 ENDIF 672 675 676 c Verification synchronize AMBIPOLAR DIFFUSION & CHEMISTRY 677 678 IF ((ok_iondiff) .and. (NINT(RDAY/dtime).ne.nbapp_chem)) THEN 679 WRITE(lunout,*)'nbapp_chem .NE. day_step' 680 WRITE(lunout,*)'nbapp_chem ', nbapp_chem 681 WRITE(lunout,*)'day_step ', NINT(RDAY/dtime) 682 WRITE(lunout,*)'nbapp_chem must be equal to day_step' 683 STOP 684 ENDIF 685 673 686 c Initialisation des sorties 674 687 c=========================== … … 1050 1063 !========= 1051 1064 1052 nbapp_chem = 240001065 1053 1066 chempas = nint(rday/pdtphys/nbapp_chem) 1054 1067 zctime = dtime*real(chempas) ! chemical timestep … … 1659 1672 1660 1673 c==================================================================== 1674 1675 c================================== 1676 ! -- ION AMBIPOLAR DIFFUSION --- 1677 c================================== 1678 1679 d_q_iondif(:,:,:)=0 1680 1681 IF (callthermos .and. ok_chem .and. ok_ionchem) THEN 1682 IF (ok_iondiff) THEN 1683 1684 call iondiff_red(klon,klev,nqmax,pplay,paprs, 1685 & t_seri,tr_seri,pphis, 1686 & gmtime,latitude_deg,longitude_deg, 1687 & pdtphys,d_t_euv,d_t_conduc,d_q_moldif, 1688 & d_q_iondif) 1689 1690 !write (*,*) 'TITI EST PASSE PAR LA' 1691 ! --- update tendencies tracers --- 1692 1693 DO iq = 1, nqmax 1694 IF (type_tr(iq) .eq. 2) THEN 1695 DO k=1,klev 1696 DO ig=1,klon 1697 tr_seri(ig,k,iq)= max(tr_seri(ig,k,iq)+ 1698 & d_q_iondif(ig,k,iq)*dtime,1.e-30) 1699 ENDDO 1700 ENDDO 1701 ENDIF 1702 ENDDO 1703 ENDIF ! ok_iondiff 1704 ENDIF ! callthermos & ok_chem & ok_ionchem 1705 1706 c==================================================================== 1661 1707 ! tests: output tendencies 1662 1708 ! call writefield_phy('physiq_dtrad',dtrad,klev) -
trunk/LMDZ.VENUS/libf/phyvenus/phytrac_chimie.F
r2836 r3035 391 391 392 392 ! solar zenith angle 393 394 sza_local = acos(cos(lat_local(ilon))*cos(lon_local(ilon)) 393 ! sza_local = acos(cos(lat_local(ilon))*cos(lon_local(ilon)) 394 ! $ *cos(lon_sun) + cos(lat_local(ilon)) 395 ! $ *sin(lon_local(ilon))*sin(lon_sun))*180./rpi 396 397 sza_local = cos(lat_local(ilon))*cos(lon_local(ilon)) 395 398 $ *cos(lon_sun) + cos(lat_local(ilon)) 396 $ *sin(lon_local(ilon))*sin(lon_sun))*180./rpi 399 $ *sin(lon_local(ilon))*sin(lon_sun) 400 401 ! Security - Handle rare cases where |sza_local| > 1 402 sza_local = min(sza_local,1.) 403 sza_local = max(-1.,sza_local) 404 sza_local = acos(sza_local)*180./rpi 397 405 398 406 ! electron temperature
Note: See TracChangeset
for help on using the changeset viewer.