Changeset 2415
- Timestamp:
- Oct 13, 2020, 5:11:50 PM (4 years ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 10 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/README
r2414 r2415 3159 3159 Correction of a bad unit for qdusttotal0 and qdusttotal1 (diagnostic outputs of the stormdust scheme) 3160 3160 + some corrections of comments about aerosol and tauref in aeropacity 3161 3162 == 13/10/2020 == EM 3163 Code cleanup: remove "tauref" and replace it with two distinct variables, 3164 "tau_pref_gcm" and "tau_pref_scenario" which are repectively the visible 3165 dust opacity column (at 610 Pa) in the GCM and prescribed by a dust scenario. -
trunk/LMDZ.MARS/libf/phymars/aeropacity_mod.F
r2414 r2415 10 10 11 11 SUBROUTINE aeropacity(ngrid,nlayer,nq,zday,pplay,pplev,ls, 12 & pq,tauscaling,tauref,tau,taucloudtes,aerosol,dsodust,reffrad, 12 & pq,tauscaling,tau_pref_scenario,tau_pref_gcm, 13 & tau,taucloudtes,aerosol,dsodust,reffrad, 13 14 & QREFvis3d,QREFir3d,omegaREFir3d, 14 15 & totstormfract,clearatm,dsords,dsotop, … … 53 54 c - reference pressure is now set to 610Pa (not 700Pa) 54 55 c 55 c input:56 c -----57 c ngrid Number of gridpoint of horizontal grid58 c nlayer Number of layer59 c nq Number of tracer60 c zday Date (time since Ls=0, in martian days)61 c ls Solar longitude (Ls) , radian62 c pplay,pplev pressure (Pa) in the middle and boundary of each layer63 c pq Dust mixing ratio (used if tracer =T and active=T).64 c reffrad(ngrid,nlayer,naerkind) Aerosol effective radius65 c QREFvis3d(ngrid,nlayer,naerkind) \ 3d extinction coefficients66 c QREFir3d(ngrid,nlayer,naerkind) / at reference wavelengths;67 c omegaREFir3d(ngrid,nlayer,naerkind) / at reference wavelengths;68 c69 c output:70 c -------71 c tauref Prescribed mean column optical depth at 610 Pa72 c tau Column total visible dust optical depth at each point73 c aerosol aerosol(ig,l,1) is the dust optical74 c depth in layer l, grid point ig75 c taualldust CW17 total opacity for all dust scatterer stormdust included76 c77 56 c======================================================================= 78 57 include "callkeys.h" … … 85 64 c Input/Output 86 65 c ------------ 87 INTEGER, INTENT(IN) :: ngrid,nlayer,nq 88 REAL, INTENT(IN) :: ls,zday 89 REAL, INTENT(IN) :: pplev(ngrid,nlayer+1),pplay(ngrid,nlayer) 90 REAL, INTENT(IN) :: pq(ngrid,nlayer,nq) 91 REAL, INTENT(OUT) :: tauref(ngrid) ! prescribed or computed 92 ! column dust opacity 93 REAL, INTENT(OUT) :: tau(ngrid,naerkind) 94 REAL,INTENT(OUT) :: taucloudtes(ngrid)! Cloud opacity at infrared 95 ! reference wavelength using 96 ! Qabs instead of Qext 97 ! (direct comparison with TES) 98 REAL, INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) 99 REAL, INTENT(OUT) :: dsodust(ngrid,nlayer) 66 INTEGER,INTENT(IN) :: ngrid ! number of atmospheric columns 67 INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers 68 INTEGER,INTENT(IN) :: nq ! number of tracers 69 REAL,INTENT(IN) :: ls ! Solar Longitude (rad) 70 REAL,INTENT(IN) :: zday ! date (in martian sols) since Ls=0 71 REAL,INTENT(IN) :: pplay(ngrid,nlayer) ! pressure (Pa) in the middle of 72 ! each atmospheric layer 73 REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! pressure (Pa) at the boundaries 74 ! of the atmospheric layers 75 REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! tracers 76 REAL,INTENT(OUT) :: tau_pref_scenario(ngrid) ! prescribed dust column 77 ! visible opacity at odpref from scenario 78 REAL,INTENT(OUT) :: tau_pref_gcm(ngrid) ! computed dust column 79 ! visible opacity at odpref in the GCM 80 REAL,INTENT(OUT) :: tau(ngrid,naerkind) ! column total visible 81 ! optical depth of each aerosol 82 REAL,INTENT(OUT) :: taucloudtes(ngrid)! Water ice cloud opacity at 83 ! infrared reference wavelength using 84 ! Qabs instead of Qext 85 ! (for direct comparison with TES) 86 REAL, INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) ! optical 87 ! depth of each aerosl in each layer 88 REAL, INTENT(OUT) :: dsodust(ngrid,nlayer) ! density scaled opacity 89 ! of (background) dust 100 90 REAL, INTENT(OUT) :: dsords(ngrid,nlayer) !dso of stormdust 101 91 REAL, INTENT(OUT) :: dsotop(ngrid,nlayer) !dso of topdust 102 REAL, INTENT(INOUT) :: reffrad(ngrid,nlayer,naerkind) 103 REAL, INTENT(IN) :: QREFvis3d(ngrid,nlayer,naerkind) 104 REAL, INTENT(IN) :: QREFir3d(ngrid,nlayer,naerkind) 105 REAL, INTENT(IN) :: omegaREFir3d(ngrid,nlayer,naerkind) 106 LOGICAL, INTENT(IN) :: clearatm 107 REAL, INTENT(IN) :: totstormfract(ngrid) 108 LOGICAL, INTENT(IN) :: nohmons 92 REAL, INTENT(INOUT) :: reffrad(ngrid,nlayer,naerkind) ! effective radius 93 ! of the aerosols in the grid boxes 94 REAL, INTENT(IN) :: QREFvis3d(ngrid,nlayer,naerkind) ! 3D extinction 95 ! coefficients (in the visible) of aerosols 96 REAL, INTENT(IN) :: QREFir3d(ngrid,nlayer,naerkind) ! 3D extinction 97 ! coefficients (in the infra-red) of aerosols 98 REAL, INTENT(IN) :: omegaREFir3d(ngrid,nlayer,naerkind) ! at the 99 ! reference wavelengths 100 LOGICAL, INTENT(IN) :: clearatm ! true to compute RT without stormdust 101 ! and false to compute RT in rocket dust storms 102 REAL, INTENT(IN) :: totstormfract(ngrid) ! mesh fraction with a rocket 103 ! dust storm 104 LOGICAL, INTENT(IN) :: nohmons ! true to compute RT without slope wind 105 ! topdust, false to compute RT in the topdust 109 106 REAL, INTENT(IN) :: alpha_hmons(ngrid) 110 107 REAL, INTENT(OUT) :: tauscaling(ngrid) ! Scaling factor for qdust and Ndust 111 REAL,INTENT(IN) :: totcloudfrac(ngrid) ! total cloud fraction 112 LOGICAL,INTENT(IN) :: clearsky ! true for part without clouds,false for part with clouds (total or sub-grid clouds) 108 REAL,INTENT(IN) :: totcloudfrac(ngrid) ! total water ice cloud fraction 109 LOGICAL,INTENT(IN) :: clearsky ! true to compute RT without water ice clouds 110 ! false to compute RT with clouds (total or sub-grid clouds) 113 111 c 114 112 c Local variables : … … 237 235 END IF ! end of if firstcall 238 236 239 ! 1. Get prescribed tauref, Dust column optical depth at "odpref" Pa 240 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 241 IF(freedust) THEN 242 tauref(:) = 0. ! tauref is computed after, instead of being forced 243 244 ELSE IF(iaervar.eq.1) THEN 237 ! 1. Get prescribed tau_pref_scenario, Dust column optical depth at "odpref" Pa 238 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 239 240 IF(iaervar.eq.1) THEN 245 241 do ig=1, ngrid 246 tau ref(ig)=max(tauvis,1.e-9) ! tauvis=cste (set in callphys.def242 tau_pref_scenario(ig)=max(tauvis,1.e-9) ! tauvis=cste (set in callphys.def 247 243 ! or read in starfi 248 244 end do 249 245 ELSE IF (iaervar.eq.2) THEN ! << "Viking" Scenario>> 250 246 251 tau ref(1) = 0.7+.3*cos(ls+80.*pi/180.) ! like seen by VL1247 tau_pref_scenario(1) = 0.7+.3*cos(ls+80.*pi/180.) ! like seen by VL1 252 248 do ig=2,ngrid 253 tau ref(ig) = tauref(1)249 tau_pref_scenario(ig) = tau_pref_scenario(1) 254 250 end do 255 251 … … 262 258 if (latitude(ig).ge.0) then 263 259 ! Northern hemisphere 264 tau ref(ig)= tauN +260 tau_pref_scenario(ig)= tauN + 265 261 & (taueq-tauN)*0.5*(1+tanh((45-latitude(ig)*180./pi)*6/60)) 266 262 else 267 263 ! Southern hemisphere 268 tau ref(ig)= tauS +264 tau_pref_scenario(ig)= tauS + 269 265 & (taueq-tauS)*0.5*(1+tanh((45+latitude(ig)*180./pi)*6/60)) 270 266 endif 271 267 enddo ! of do ig=1,ngrid 272 268 ELSE IF (iaervar.eq.5) THEN ! << Escalier Scenario>> 273 tau ref(1) = 2.5269 tau_pref_scenario(1) = 2.5 274 270 if ((ls.ge.30.*pi/180.).and.(ls.le.150.*pi/180.)) 275 & tau ref(1) = .2271 & tau_pref_scenario(1) = .2 276 272 277 273 do ig=2,ngrid 278 tau ref(ig) = tauref(1)274 tau_pref_scenario(ig) = tau_pref_scenario(1) 279 275 end do 280 276 ELSE IF ((iaervar.ge.6).and.(iaervar.le.8)) THEN 281 277 ! clim, cold or warm synthetic scenarios 282 call read_dust_scenario(ngrid,nlayer,zday,pplev,tauref) 278 call read_dust_scenario(ngrid,nlayer,zday,pplev, 279 & tau_pref_scenario) 283 280 ELSE IF ((iaervar.ge.24).and.(iaervar.le.34)) 284 281 & THEN ! << MY... dust scenarios >> 285 call read_dust_scenario(ngrid,nlayer,zday,pplev,tauref) 282 call read_dust_scenario(ngrid,nlayer,zday,pplev, 283 & tau_pref_scenario) 286 284 ELSE IF ((iaervar.eq.4).or. 287 285 & ((iaervar.ge.124).and.(iaervar.le.126))) THEN 288 286 ! "old" TES assimation dust scenario (values at 700Pa in files!) 289 call read_dust_scenario(ngrid,nlayer,zday,pplev,tauref) 287 call read_dust_scenario(ngrid,nlayer,zday,pplev, 288 & tau_pref_scenario) 290 289 ELSE 291 290 call abort_physic("aeropacity","wrong value for iaervar",1) … … 354 353 DO ig=1,ngrid 355 354 zp=odpref/pplay(ig,l) 356 aerosol(ig,l,1)= tau ref(ig)/odpref *355 aerosol(ig,l,1)= tau_pref_scenario(ig)/odpref * 357 356 s (pplev(ig,l)-pplev(ig,l+1)) 358 357 s *max( exp(.03*(1.-max(zp,1.))) , 1.E-3 ) … … 589 588 DO l=1,nlayer 590 589 DO ig=1,ngrid 591 tau ref(ig) = tauref(ig) +590 tau_pref_gcm(ig) = tau_pref_gcm(ig) + 592 591 & aerosol(ig,l,iaerdust(iaer)) 593 592 ENDDO 594 593 ENDDO 595 594 ENDDO 596 tau ref(:) = tauref(:) * odpref / pplev(:,1)595 tau_pref_gcm(:) = tau_pref_gcm(:) * odpref / pplev(:,1) 597 596 598 597 c-------------------------------------------------- … … 656 655 !! -- the additional reference opacity will 657 656 !! thus be taulocref*odpref/pplev 658 tauuser(ig)=max( tau ref(ig) * pplev(ig,1) /odpref ,657 tauuser(ig)=max( tau_pref_gcm(ig) * pplev(ig,1) /odpref , 659 658 & taulocref * yeah ) 660 659 … … 696 695 !! Mass mixing ratio perturbation due to local dust storm in each layer 697 696 more_dust(ig,l,1)= 698 & (tauuser(ig)-(tau ref(ig)697 & (tauuser(ig)-(tau_pref_gcm(ig) 699 698 & * pplev(ig,1) /odpref)) / 700 699 & int_factor(ig) 701 700 more_dust(ig,l,2)= 702 & (tauuser(ig)-(tau ref(ig) *701 & (tauuser(ig)-(tau_pref_gcm(ig) * 703 702 & pplev(ig,1) /odpref)) 704 703 & / int_factor(ig) * … … 716 715 ENDDO 717 716 ENDIF !! IF (localstorm) 718 tau ref(:)=0.717 tau_pref_gcm(:)=0. 719 718 ENDIF !! IF (firstcall) 720 719 #endif … … 723 722 ! 3.1. Compute "tauscaling", the dust rescaling coefficient and adjust 724 723 ! aerosol() dust opacities accordingly 725 call compute_tauscaling(ngrid,nlayer,naerkind,naerdust, 726 & zday,pplev,tauref,tauscaling,aerosol)727 728 ! 3.2. Recompute tau ref, the reference dust opacity, based on dust tracer724 call compute_tauscaling(ngrid,nlayer,naerkind,naerdust,zday,pplev, 725 & tau_pref_scenario,tauscaling,aerosol) 726 727 ! 3.2. Recompute tau_pref_gcm, the reference dust opacity, based on dust tracer 729 728 ! mixing ratios and their optical properties 730 729 731 730 IF (freedust) THEN 732 ! tauref has been initialized to 0 before. 731 ! Initialisation : 732 tau_pref_gcm(:)=0 733 733 DO iaer=1,naerdust 734 734 DO l=1,nlayer … … 744 744 ENDIF 745 745 #endif 746 c MV19: tau refmust ALWAYS contain the opacity of all dust tracers747 ! GCM DUST OPTICAL DEPTH tau refis to be compared748 ! with the observation CDOD tau ref_scenario746 c MV19: tau_pref_gcm must ALWAYS contain the opacity of all dust tracers 747 ! GCM DUST OPTICAL DEPTH tau_pref_gcm is to be compared 748 ! with the observation CDOD tau_pref_scenario 749 749 ! => visible wavelength 750 750 IF (name_iaer(iaerdust(iaer)).eq."dust_doubleq") THEN 751 tau ref(ig) = tauref(ig) +751 tau_pref_gcm(ig) = tau_pref_gcm(ig) + 752 752 & ( 0.75 * QREFvis3d(ig,l,iaerdust(iaer)) / 753 753 & ( rho_dust * reffrad(ig,l,iaerdust(iaer)) ) ) * … … 755 755 & ( pplev(ig,l) - pplev(ig,l+1) ) / g 756 756 ELSE IF (name_iaer(iaerdust(iaer)).eq."stormdust_doubleq") THEN 757 tau ref(ig) = tauref(ig) +757 tau_pref_gcm(ig) = tau_pref_gcm(ig) + 758 758 & ( 0.75 * QREFvis3d(ig,l,iaerdust(iaer)) / 759 759 & ( rho_dust * reffrad(ig,l,iaerdust(iaer)) ) ) * … … 761 761 & ( pplev(ig,l) - pplev(ig,l+1) ) / g 762 762 ELSE IF (name_iaer(iaerdust(iaer)).eq."topdust_doubleq") THEN 763 tau ref(ig) = tauref(ig) +763 tau_pref_gcm(ig) = tau_pref_gcm(ig) + 764 764 & ( 0.75 * QREFvis3d(ig,l,iaerdust(iaer)) / 765 765 & ( rho_dust * reffrad(ig,l,iaerdust(iaer)) ) ) * … … 771 771 ENDDO 772 772 ENDDO 773 tauref(:) = tauref(:) * odpref / pplev(:,1) 773 tau_pref_gcm(:) = tau_pref_gcm(:) * odpref / pplev(:,1) 774 ELSE 775 ! dust opacity strictly follows what is imposed by the dust scenario 776 tau_pref_gcm(:)=tau_pref_scenario(:) 774 777 ENDIF ! of IF (freedust) 775 778 -
trunk/LMDZ.MARS/libf/phymars/callradite_mod.F
r2409 r2415 8 8 $ emis,mu0,pplev,pplay,pt,tsurf,fract,dist_sol,igout, 9 9 $ dtlw,dtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw, 10 $ fluxtop_sw,tauref,tau,aerosol,dsodust,tauscaling, 10 $ fluxtop_sw,tau_pref_scenario,tau_pref_gcm, 11 & tau,aerosol,dsodust,tauscaling, 11 12 $ taucloudtes,rdust,rice,nuice,co2ice,rstormdust,rtopdust, 12 13 $ totstormfract,clearatm,dsords,dsotop,alpha_hmons,nohmons, … … 153 154 c fluxtop_sw(ngrid,2) outgoing upward flux SW for solar band#2 (W.m-2) 154 155 155 c tauref Prescribed mean column optical depth at 610 Pa156 156 c tau Column total visible dust optical depth at each point 157 157 c aerosol(ngrid,nlayer,naerkind) aerosol extinction optical depth … … 187 187 REAL,INTENT(OUT) :: fluxsurf_lw(ngrid), fluxtop_lw(ngrid) 188 188 REAL,INTENT(OUT) :: fluxsurf_sw(ngrid,2), fluxtop_sw(ngrid,2) 189 190 REAL,INTENT(OUT) :: tauref(ngrid), tau(ngrid,naerkind) 189 REAL,INTENT(OUT) :: tau_pref_scenario(ngrid) ! prescribed dust column 190 ! visible opacity at odpref from scenario 191 REAL,INTENT(OUT) :: tau_pref_gcm(ngrid) ! computed dust column 192 ! visible opacity at odpref in the GCM 193 REAL,INTENT(OUT) :: tau(ngrid,naerkind) 191 194 REAL,INTENT(OUT) :: taucloudtes(ngrid)! Cloud opacity at infrared 192 195 ! reference wavelength using … … 204 207 REAL,INTENT(IN) :: totstormfract(ngrid) ! dust storm mesh fraction 205 208 REAL,INTENT(OUT) :: rstormdust(ngrid,nlayer) ! Storm dust geometric mean radius (m) 206 REAL,INTENT( INOUT) :: dsords(ngrid,nlayer) ! density scaled opacity for rocket dust storm dust209 REAL,INTENT(OUT) :: dsords(ngrid,nlayer) ! density scaled opacity for rocket dust storm dust 207 210 208 211 c entrainment by slope wind … … 210 213 REAL, INTENT(IN) :: alpha_hmons(ngrid) ! sub-grid scale topography mesh fraction 211 214 REAL,INTENT(OUT) :: rtopdust(ngrid,nlayer) ! Topdust geometric mean radius (m) 212 REAL,INTENT( INOUT) :: dsotop(ngrid,nlayer) ! density scaled opacity for topmons dust215 REAL,INTENT(OUT) :: dsotop(ngrid,nlayer) ! density scaled opacity for topmons dust 213 216 214 217 c sub-grid scale water ice clouds … … 423 426 c Computing aerosol optical depth in each layer: 424 427 CALL aeropacity(ngrid,nlayer,nq,zday,pplay,pplev,ls, 425 & pq,tauscaling,tauref,tau,taucloudtes,aerosol,dsodust,reffrad, 428 & pq,tauscaling,tau_pref_scenario,tau_pref_gcm, 429 & tau,taucloudtes,aerosol,dsodust,reffrad, 426 430 & QREFvis3d,QREFir3d,omegaREFir3d, 427 431 & totstormfract,clearatm,dsords,dsotop, -
trunk/LMDZ.MARS/libf/phymars/compute_dtau_mod.F90
r2409 r2415 5 5 REAL,SAVE :: ti_injection_sol ! time of beginning injection 6 6 REAL,SAVE :: tf_injection_sol ! time of end injection 7 REAL,PARAMETER :: t_scenario_sol=14/24. ! time of tauref_scenario 7 REAL,PARAMETER :: t_scenario_sol=14/24. ! time of day at which 8 ! tau_pref_scenario is deemed exact 8 9 9 10 REAL,SAVE,ALLOCATABLE :: dtau(:) ! Dust opacity difference (at 610Pa) … … 13 14 14 15 SUBROUTINE compute_dtau(ngrid,nlayer, & 15 zday,pplev,tau ref,&16 zday,pplev,tau_pref_gcm, & 16 17 ptimestep,dustliftday,local_time) 17 18 … … 31 32 REAL, INTENT(in) :: zday ! date at lon=0, in fraction of sols 32 33 REAL, INTENT(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa) 33 REAL, INTENT(in) :: tauref(ngrid) ! Computed dust opacity at 610Pa 34 REAL, INTENT(in) :: tau_pref_gcm(ngrid) ! Visible dust opacity column 35 ! at 610Pa as computed in the GCM 34 36 REAL, INTENT(in) :: ptimestep 35 37 REAL, INTENT(in) :: local_time(ngrid) … … 38 40 INTEGER :: ig, l 39 41 INTEGER, SAVE :: nb_daystep ! nomber of step a day 40 REAL :: tauref_scenario(ngrid) ! from dust scenario 42 REAL :: tau_pref_target(ngrid) ! dust opacity column at odpref=610 Pa 43 ! as extracted from dust scenario 41 44 REAL :: zday_scenario 42 45 REAL,ALLOCATABLE,SAVE :: local_time_prev(:) … … 59 62 ENDIF 60 63 61 ! 1. Obtain tau ref_scenariofrom dust scenario at zday+164 ! 1. Obtain tau_pref_target from dust scenario at zday+1 62 65 if (iaervar.eq.1) then 63 tau ref_scenario= tauvis66 tau_pref_target = tauvis 64 67 else 65 68 zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario opacity is measured at 14:00 66 69 zday_scenario=zday_scenario+1 ! opacity of the dust scenario is read the day after 67 70 call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev, & 68 tau ref_scenario)71 tau_pref_target) 69 72 endif 70 73 ! for diagnostics 71 call WRITEDIAGFI(ngrid,"tauref_scenario","tauref_scenario", & 72 "",2,tauref_scenario) 74 call WRITEDIAGFI(ngrid,"tau_pref_target", & 75 "target visible dust opacity column at 610Pa", & 76 "",2,tau_pref_target) 73 77 74 78 ! 2. Compute dtau() and dustliftday() … … 76 80 IF ((local_time(ig).ge.t_scenario_sol).and. & 77 81 (local_time_prev(ig).lt.(t_scenario_sol)))THEN 78 dtau(ig)=tau ref_scenario(ig)-tauref(ig)82 dtau(ig)=tau_pref_target(ig)-tau_pref_gcm(ig) 79 83 ENDIF 80 84 81 ! Use dtau (when positi ove) to compute dustliftday85 ! Use dtau (when positive) to compute dustliftday 82 86 IF (dtau(ig).LT.0) THEN 83 87 dustliftday(ig)=0. -
trunk/LMDZ.MARS/libf/phymars/dust_param_mod.F90
r2413 r2415 12 12 13 13 REAL,PARAMETER :: odpref = 610. ! Reference pressure (Pa) of 14 ! DOD (Dust optical Depth) tau ref14 ! DOD (Dust optical Depth) tau_pref_* 15 15 16 16 REAL,SAVE,ALLOCATABLE :: tauscaling(:) ! Convertion factor for qdust and Ndust -
trunk/LMDZ.MARS/libf/phymars/dust_scaling_mod.F90
r2413 r2415 7 7 subroutine compute_tauscaling(ngrid,nlayer,naerkind,naerdust, & 8 8 zday,pplev, & 9 tau ref,tauscaling,aerosol)9 tau_pref_scenario,tauscaling,aerosol) 10 10 11 11 use dust_param_mod, only: tauscaling_mode, odpref … … 20 20 real,intent(in) :: zday 21 21 real,intent(in) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa) 22 real,intent(in) :: tauref(ngrid) ! prescribed dust opacity 22 real,intent(in) :: tau_pref_scenario(ngrid) ! prescribed visible dust 23 ! opacity column at odpref reference pressure 23 24 real,intent(out) :: tauscaling(ngrid) ! dust scaling factor 24 25 real,intent(inout) :: aerosol(ngrid,nlayer,naerkind) ! opacities … … 36 37 if (tauscaling_mode == 1) then 37 38 ! GCM v5.3 style: tauscaling is computed so that 38 ! aerosol() opacities correspond to the prescribed tau ref()39 ! aerosol() opacities correspond to the prescribed tau_pref_scenario() 39 40 40 41 ! 1. compute dust column opacity using aerosol() dusts … … 49 50 50 51 ! 2. compute the scaling factor 51 tauscaling(:)=tau ref(:)*pplev(:,1)/odpref/taudust(:)52 tauscaling(:)=tau_pref_scenario(:)*pplev(:,1)/odpref/taudust(:) 52 53 endif ! of if (tauscaling_mode == 1) 53 54 -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r2414 r2415 281 281 REAL fluxtop_lw(ngrid) !Outgoing LW (IR) flux to space (W.m-2) 282 282 REAL fluxtop_sw(ngrid,2) !Outgoing SW (solar) flux to space (W.m-2) 283 REAL tauref(ngrid) ! Reference column optical depth at odpref 283 REAL tau_pref_scenario(ngrid) ! prescribed dust column visible opacity 284 ! at odpref 285 REAL tau_pref_gcm(ngrid) ! dust column visible opacity at odpref in the GCM 284 286 c rocket dust storm 285 287 REAL totstormfract(ngrid) ! fraction of the mesh where the dust storm is contained … … 928 930 & emis,mu0,zplev,zplay,pt,tsurf,fract,dist_sol,igout, 929 931 & zdtlw,zdtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw, 930 & fluxtop_sw,tauref,tau,aerosol,dsodust,tauscaling, 932 & fluxtop_sw,tau_pref_scenario,tau_pref_gcm, 933 & tau,aerosol,dsodust,tauscaling, 931 934 & taucloudtes,rdust,rice,nuice,co2ice,rstormdust,rtopdust, 932 935 & totstormfract,clearatm,dsords,dsotop,alpha_hmons,nohmons, … … 943 946 & albedo,emis,mu0,zplev,zplay,pt,tsurf,fract, 944 947 & dist_sol,igout,zdtlwclf,zdtswclf,fluxsurf_lwclf, 945 & fluxsurf_swclf,fluxtop_lwclf,fluxtop_swclf,tauref, 948 & fluxsurf_swclf,fluxtop_lwclf,fluxtop_swclf, 949 & tau_pref_scenario,tau_pref_gcm, 946 950 & tau,aerosol,dsodust,tauscaling,taucloudtesclf,rdust, 947 951 & rice,nuice,co2ice,rstormdust,rtopdust,totstormfract, … … 999 1003 & longitude(igout)*180/pi 1000 1004 1001 write(*,'(" tau ref(",f4.0," Pa) =",f9.6,1005 write(*,'(" tau_pref_gcm(",f4.0," Pa) =",f9.6, 1002 1006 & " tau(",f4.0," Pa) =",f9.6)') 1003 & odpref,tau ref(igout),1007 & odpref,tau_pref_gcm(igout), 1004 1008 & odpref,tau(igout,1)*odpref/zplev(igout,1) 1005 1009 c --------------------------------------------------------- … … 1140 1144 c output 1141 1145 & pdqrds,wspeed,dsodust,dsords,dsotop, 1142 & tau ref)1146 & tau_pref_scenario,tau_pref_gcm) 1143 1147 1144 1148 c update the tendencies of both dust after vertical transport … … 1204 1208 & nohmons,hsummit, 1205 1209 & pdqtop,wtop,dsodust,dsords,dsotop, 1206 & tau ref)1210 & tau_pref_scenario,tau_pref_gcm) 1207 1211 1208 1212 … … 1231 1235 1232 1236 CALL compute_dtau(ngrid,nlayer, 1233 & zday,pplev,tau ref,1237 & zday,pplev,tau_pref_gcm, 1234 1238 & ptimestep,dustliftday,local_time) 1235 1239 endif ! end of if (dustinjection.gt.0) … … 2495 2499 call wstats(ngrid,"watercap","H2O ice cover", 2496 2500 & "kg.m-2",2,watercap) 2497 call wstats(ngrid,"tauref","reference dod at 610 Pa","NU", 2498 & 2,tauref) 2501 call wstats(ngrid,"tau_pref_scenario", 2502 & "prescribed visible dod at 610 Pa","NU", 2503 & 2,tau_pref_scenario) 2504 call wstats(ngrid,"tau_pref_gcm", 2505 & "visible dod at 610 Pa in the GCM","NU", 2506 & 2,tau_pref_gcm) 2499 2507 call wstats(ngrid,"fluxsurf_lw", 2500 2508 & "Thermal IR radiative flux to surface","W.m-2",2, … … 2752 2760 comm_SWDOWNZ(1:ngrid) = fluxsurf_sw_tot(1:ngrid) 2753 2761 !state real TAU_DUST ij misc 1 - h "TAU_DUST" "REFERENCE VISIBLE DUST OPACITY" "" 2754 comm_TAU_DUST(1:ngrid) = tau ref(1:ngrid)2762 comm_TAU_DUST(1:ngrid) = tau_pref_gcm(1:ngrid) 2755 2763 !state real RDUST ikj misc 1 - h "RDUST" "DUST RADIUS" "m" 2756 2764 comm_RDUST(1:ngrid,1:nlayer) = rdust(1:ngrid,1:nlayer) … … 3086 3094 c ---------------------------------------------------------- 3087 3095 3088 call WRITEDIAGFI(ngrid,'tauref', 3089 & 'Dust ref opt depth','NU',2,tauref) 3096 call WRITEDIAGFI(ngrid,'tau_pref_scenario', 3097 & 'Prescribed visible dust optical depth at 610Pa', 3098 & 'NU',2,tau_pref_scenario) 3099 3100 call WRITEDIAGFI(ngrid,'tau_pref_gcm', 3101 & 'Visible dust optical depth at 610Pa in the GCM', 3102 & 'NU',2,tau_pref_gcm) 3090 3103 3091 3104 if (tracer.and.(dustbin.ne.0)) then … … 3382 3395 3383 3396 #ifndef MESOSCALE 3384 write(*,'("Ls =",f11.6," tauref(",f4.0," Pa) =",f9.6)') 3385 & zls*180./pi,odpref,tauref 3397 write(*, 3398 & '("Ls =",f11.6," tau_pref_scenario(",f4.0," Pa) =",f9.6)') 3399 & zls*180./pi,odpref,tau_pref_scenario 3386 3400 c ---------------------------------------------------------------------- 3387 3401 c Output in grads file "g1d" (ONLY when using testphys1d) … … 3491 3505 & 'total mass of dust only', 3492 3506 & 'kg.m-2',0,mdusttot) 3493 call WRITEDIAGFI(ngrid,'tauref', 3494 & 'Dust ref opt depth','NU',0,tauref) 3507 call WRITEDIAGFI(ngrid,'tau_pref_scenario', 3508 & 'Prescribed dust ref opt depth at 610 Pa', 3509 & 'NU',0,tau_pref_scenario) 3510 call WRITEDIAGFI(ngrid,'tau_pref_gcm', 3511 & 'Dust ref opt depth at 610 Pa in the GCM', 3512 & 'NU',0,tau_pref_gcm) 3495 3513 call WRITEDIAGFI(ngrid,'rdsdqsdust', 3496 3514 & 'deposited surface stormdust mass', -
trunk/LMDZ.MARS/libf/phymars/read_dust_scenario.F90
r2398 r2415 1 subroutine read_dust_scenario(ngrid,nlayer,zday,pplev,tau ref)1 subroutine read_dust_scenario(ngrid,nlayer,zday,pplev,tau_pref_scenario) 2 2 3 3 ! Reading of the dust scenario file … … 6 6 use geometry_mod, only: latitude, longitude ! in radians 7 7 use datafile_mod, only: datadir 8 use dust_param_mod, only: odpref 8 9 implicit none 9 10 … … 13 14 real, intent(in) :: zday ! date in martian days 14 15 real, dimension(ngrid,nlayer+1), intent(in) :: pplev 15 real, dimension(ngrid), intent(out) :: tauref 16 real, dimension(ngrid), intent(out) :: tau_pref_scenario ! visible dust column 17 ! opacity at odpref from scenario 16 18 17 19 ! Local variables … … 333 335 tau= dlat*(dlon*(tau1(2)+tau1(3)-tau1(1)-tau1(4))+tau1(1)-tau1(3)) +dlon*(tau1(4)-tau1(3))+tau1(3) 334 336 335 tau ref(ig)=tau337 tau_pref_scenario(ig)=tau 336 338 ! 337 339 enddo ! of ig=1,ngrid … … 342 344 ! - the 1.3 conversion factor from IR absorbtion opacity to 343 345 ! IR extinction opacity 344 tau ref(:)=tauref(:)*1.3*(610./700.)346 tau_pref_scenario(:)=tau_pref_scenario(:)*1.3*(odpref/700.) 345 347 endif 346 348 … … 349 351 ! needed to decrease opacity (*0.825) to compensate overestimation of 350 352 ! heating rates 351 tau ref(:)=tauref(:)*0.825/1.3353 tau_pref_scenario(:)=tau_pref_scenario(:)*0.825/1.3 352 354 endif 353 355 -
trunk/LMDZ.MARS/libf/phymars/rocketduststorm_mod.F90
r2252 r2415 31 31 ! output 32 32 pdqrds,wrad,dsodust,dsords,dsotop, & 33 tau ref)33 tau_pref_scenario,tau_pref_gcm) 34 34 35 35 USE tracer_mod, only: igcm_stormdust_mass,igcm_stormdust_number & … … 92 92 REAL, INTENT(OUT) :: pdqrds(ngrid,nlayer,nq) ! tendancy field for dust when detraining 93 93 REAL, INTENT(OUT) :: wrad(ngrid,nlayer+1) ! vertical speed within the rocket dust storm 94 REAL, INTENT(INOUT) :: dsodust(ngrid,nlayer) ! density scaled opacity of env. dust 95 REAL, INTENT(INOUT) :: dsords(ngrid,nlayer) ! density scaled opacity of storm dust 96 REAL, INTENT(INOUT) :: dsotop(ngrid,nlayer) ! density scaled opacity of topmons dust 97 REAL, INTENT(OUT) :: tauref(ngrid) 98 94 REAL, INTENT(OUT) :: dsodust(ngrid,nlayer) ! density scaled opacity of env. dust 95 REAL, INTENT(OUT) :: dsords(ngrid,nlayer) ! density scaled opacity of storm dust 96 REAL, INTENT(OUT) :: dsotop(ngrid,nlayer) ! density scaled opacity of topmons dust 97 REAL,INTENT(OUT) :: tau_pref_scenario(ngrid) ! prescribed dust column 98 ! visible opacity at odpref 99 REAL,INTENT(OUT) :: tau_pref_gcm(ngrid) ! dust column visible opacity at 100 ! odpref in the GCM 99 101 !-------------------------------------------------------- 100 102 ! Local variables … … 246 248 emis,mu0,pplev,pplay,pt,tsurf,fract,dist_sol,igout, & 247 249 zdtlw1,zdtsw1,fluxsurf_lw1,fluxsurf_sw1,fluxtop_lw1, & 248 fluxtop_sw1,tauref,tau,aerosol,dsodust,tauscaling, & 250 fluxtop_sw1,tau_pref_scenario,tau_pref_gcm, & 251 tau,aerosol,dsodust,tauscaling, & 249 252 taucloudtes,rdust,rice,nuice,co2ice,rstormdust,rtopdust, & 250 253 totstormfract,clearatm,dsords,dsotop,alpha_hmons,nohmons,& -
trunk/LMDZ.MARS/libf/phymars/topmons_mod.F90
r2255 r2415 33 33 ! output 34 34 pdqtop,wfin,dsodust,dsords,dsotop, & 35 tau ref)35 tau_pref_scenario,tau_pref_gcm) 36 36 37 37 USE tracer_mod, only: igcm_topdust_mass,igcm_topdust_number & … … 96 96 ! output for second radiative transfer 97 97 REAL, INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) 98 REAL, INTENT(INOUT) :: dsodust(ngrid,nlayer) 99 REAL, INTENT(INOUT) :: dsords(ngrid,nlayer) 100 REAL, INTENT(INOUT) :: dsotop(ngrid,nlayer) 101 REAL, INTENT(OUT) :: tauref(ngrid) 98 REAL, INTENT(OUT) :: dsodust(ngrid,nlayer) 99 REAL, INTENT(OUT) :: dsords(ngrid,nlayer) 100 REAL, INTENT(OUT) :: dsotop(ngrid,nlayer) 101 REAL,INTENT(OUT) :: tau_pref_scenario(ngrid) ! prescribed dust column 102 ! visible opacity at odpref 103 REAL,INTENT(OUT) :: tau_pref_gcm(ngrid) ! dust column visible opacity at 104 ! odpref in the GCM 102 105 103 106 !-------------------------------------------------------- … … 271 274 emis,mu0,pplev,pplay,pt,tsurf,fract,dist_sol,igout, & 272 275 zdtlw1,zdtsw1,fluxsurf_lw1,fluxsurf_sw1,fluxtop_lw1, & 273 fluxtop_sw1,tauref,tau,aerosol,dsodust,tauscaling, & 276 fluxtop_sw1,tau_pref_scenario,tau_pref_gcm, & 277 tau,aerosol,dsodust,tauscaling, & 274 278 taucloudtes,rdust,rice,nuice,co2ice,rstormdust,rtopdust, & 275 279 totstormfract,clearatm,dsords,dsotop,alpha_hmons,nohmons,&
Note: See TracChangeset
for help on using the changeset viewer.