Changeset 3901
- Timestamp:
- Aug 20, 2025, 4:25:12 PM (2 days ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 1 deleted
- 11 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/changelog.txt
r3886 r3901 4945 4945 Fix newstart, now that initracer is a module. 4946 4946 4947 == 20/08/2025 == EM 4948 Some code tidying: 4949 - turn aeroptproperties, albedocaps, cvmgp and convadj into modules 4950 - remove useless check in callradite 4951 - clean module vlz_fi (remove obsolete #ifdef CRAY cpp alternatives) 4952 - remove function cvmgp since it was only called under #ifdef CRAY in vlz_fi -
trunk/LMDZ.MARS/libf/phymars/aeroptproperties.F
r3726 r3901 1 MODULE aeroptproperties_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 SUBROUTINE aeroptproperties(ngrid,nlayer,reffrad,nueffrad, 2 8 & QVISsQREF3d,omegaVIS3d,gVIS3d, … … 1345 1351 c================================================================== 1346 1352 1347 RETURN 1348 END 1353 END SUBROUTINE aeroptproperties 1354 1355 END MODULE aeroptproperties_mod -
trunk/LMDZ.MARS/libf/phymars/albedocaps.F90
r3726 r3901 1 module albedocaps_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine albedocaps(zls,ngrid,piceco2,piceco2_peren,psolaralb,emisref) 2 8 … … 518 524 end subroutine TES_icecap_albedo 519 525 526 end module albedocaps_mod -
trunk/LMDZ.MARS/libf/phymars/callradite_mod.F
r3756 r3901 17 17 use aeropacity_mod, only: aeropacity 18 18 use updatereffrad_mod, only: updatereffrad 19 use aeroptproperties_mod, only: aeroptproperties 19 20 use dimradmars_mod, only: ndomainsz, nflev, nsun, nir 20 21 use dimradmars_mod, only: naerkind, name_iaer, … … 288 289 c --------------------- 289 290 290 real zco2 ! volume fraction of CO2 in Mars atmosphere 291 !$OMP THREADPRIVATE(zco2) 292 DATA zco2/0.95/ 293 SAVE zco2 294 295 LOGICAL firstcall 291 LOGICAL,SAVE :: firstcall = .true. 296 292 !$OMP THREADPRIVATE(firstcall) 297 DATA firstcall/.true./298 SAVE firstcall299 293 300 294 … … 413 407 414 408 ! CALL SULW ! this step is now done in ini_yomlw_h 415 416 if (ngrid .EQ. 1) then417 if (ndomainsz .NE. 1) then418 print*419 print*,'ATTENTION !!!'420 print*,'pour tourner en 1D, '421 print*,'fixer ndomainsz=1 dans phymars/dimradmars_h'422 print*423 call exit(1)424 endif425 endif426 409 427 410 firstcall=.false. -
trunk/LMDZ.MARS/libf/phymars/co2condens_mod.F
r3739 r3901 40 40 use callkeys_mod, only: caps, co2clouds 41 41 use co2snow_mod, only: co2snow 42 use albedocaps_mod, only: albedocaps 42 43 43 44 IMPLICIT NONE -
trunk/LMDZ.MARS/libf/phymars/convadj.F
r3726 r3901 1 module convadj_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine convadj(ngrid,nlay,nq,ptimestep, 2 8 & pplay,pplev,ppopsk,lmax_th, … … 381 387 & -zq(1:ngrid,1:nlay,1:nq))/ptimestep 382 388 383 end 389 end subroutine convadj 384 390 391 end module convadj_mod -
trunk/LMDZ.MARS/libf/phymars/cvmgt.F
r38 r3901 1 FUNCTION cvmgt(x1,x2,l) 1 MODULE cvmgt_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 7 REAL FUNCTION cvmgt(x1,x2,l) 8 ! function which returns x1 if l == .true. , x2 otherwise 2 9 IMPLICIT NONE 3 10 4 REAL x1,x2,cvmgt5 LOGICAL l11 REAL,INTENT(IN) :: x1,x2 12 LOGICAL,INTENT(IN) :: l 6 13 7 14 IF(l) then … … 11 18 ENDIF 12 19 13 RETURN 14 END 20 END FUNCTION cvmgt 15 21 C 22 END MODULE cvmgt_mod -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r3855 r3901 20 20 use co2cloud_mod, only: co2cloud 21 21 use callradite_mod, only: callradite 22 use convadj_mod, only: convadj 22 23 use callsedim_mod, only: callsedim 23 24 use rocketduststorm_mod, only: rocketduststorm, dustliftday -
trunk/LMDZ.MARS/libf/phymars/swr_fouquart.F
r3759 r3901 15 15 use callkeys_mod, only: rayleigh 16 16 use swrayleigh_mod, only: swrayleigh 17 use cvmgt_mod, only: cvmgt 17 18 18 19 IMPLICIT NONE … … 112 113 S , ZTRA1(NDLON,NFLEV+1), ZTRA2(NDLON,NFLEV+1) 113 114 114 c Function115 c --------116 real CVMGT117 118 115 C -------------------------------- 119 116 C OPTICAL PARAMETERS FOR AEROSOLS … … 138 135 ZTAUAZ(JL,JK) = 0. 139 136 END DO 140 DO 106JAE=1,naerkind141 DO 105JL = 1 , KDLON137 DO JAE=1,naerkind 138 DO JL = 1 , KDLON 142 139 c Mean Extinction optical depth in the spectral band 143 140 c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ … … 154 151 S QVISsQREF3d(JL,JK,KNU,JAE)* 155 152 & omegaVIS3d(JL,JK,KNU,JAE)*gVIS3d(JL,JK,KNU,JAE) 156 105 CONTINUE157 106 CONTINUE153 ENDDO 154 ENDDO 158 155 END DO 159 156 C 160 157 DO JK = 1 , nlaylte 161 158 DO JL = 1 , KDLON 159 ! NB: function CVMGT(x1,x2,l) returns x1 if l==.true. 160 ! or x2 otherwise. Maybe worth inlining someday? 162 161 ZCGAZ(JL,JK) = CVMGT( 0., ZCGAZ(JL,JK) / ZPIZAZ(JL,JK), 163 162 S (ZPIZAZ(JL,JK).EQ.0) ) -
trunk/LMDZ.MARS/libf/phymars/swr_toon.F
r3740 r3901 15 15 use callkeys_mod, only: rayleigh 16 16 use swrayleigh_mod, only: swrayleigh 17 use cvmgt_mod, only: cvmgt 17 18 18 19 IMPLICIT NONE … … 75 76 C ARGUMENTS 76 77 C --------- 77 INTEGER KDLON, KFLEV, KNU78 REAL aerosol(NDLO2,KFLEV,naerkind), albedo(NDLO2,2),79 SPDSIG(NDLO2,KFLEV),PSEC(NDLO2)80 81 REAL QVISsQREF3d(NDLO2,KFLEV,nsun,naerkind)82 REAL omegaVIS3d(NDLO2,KFLEV,nsun,naerkind)83 REAL gVIS3d(NDLO2,KFLEV,nsun,naerkind)84 85 REAL PPSOL(NDLO2)86 REAL PFD(NDLO2,KFLEV+1),PFU(NDLO2,KFLEV+1)87 REAL PRMU(NDLO2)78 INTEGER,INTENT(IN) :: KDLON, KFLEV, KNU 79 REAL,INTENT(IN) :: aerosol(NDLO2,KFLEV,naerkind), albedo(NDLO2,2) 80 REAL,INTENT(IN) :: PDSIG(NDLO2,KFLEV),PSEC(NDLO2) 81 82 REAL,INTENT(IN) :: QVISsQREF3d(NDLO2,KFLEV,nsun,naerkind) 83 REAL,INTENT(IN) :: omegaVIS3d(NDLO2,KFLEV,nsun,naerkind) 84 REAL,INTENT(IN) :: gVIS3d(NDLO2,KFLEV,nsun,naerkind) 85 86 REAL,INTENT(IN) :: PPSOL(NDLO2) 87 REAL,INTENT(OUT) :: PFD(NDLO2,KFLEV+1),PFU(NDLO2,KFLEV+1) 88 REAL,INTENT(IN) :: PRMU(NDLO2) 88 89 89 90 C LOCAL ARRAYS … … 104 105 c End part added by Tran The Trung 105 106 106 c Function107 c --------108 real CVMGT109 110 107 c Computing TOTAL single scattering parameters by adding 111 108 c properties of all the NAERKIND kind of aerosols … … 117 114 ZTAUAZ(JL,JK) = 0. 118 115 END DO 119 DO 106JAE=1,naerkind120 DO 105JL = 1 , KDLON116 DO JAE=1,naerkind 117 DO JL = 1 , KDLON 121 118 c Mean Extinction optical depth in the spectral band 122 119 c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ … … 133 130 S QVISsQREF3d(JL,JK,KNU,JAE)* 134 131 & omegaVIS3d(JL,JK,KNU,JAE)*gVIS3d(JL,JK,KNU,JAE) 135 105 CONTINUE136 106 CONTINUE132 ENDDO 133 ENDDO 137 134 END DO 138 135 C 139 136 DO JK = 1 , nlaylte 140 137 DO JL = 1 , KDLON 138 ! NB: function CVMGT(x1,x2,l) returns x1 if l==.true. 139 ! or x2 otherwise. Maybe worth inlining someday? 141 140 ZCGAZ(JL,JK) = CVMGT( 0., ZCGAZ(JL,JK) / ZPIZAZ(JL,JK), 142 141 S (ZPIZAZ(JL,JK).EQ.0) ) -
trunk/LMDZ.MARS/libf/phymars/vlz_fi.F
r3727 r3901 10 10 c 11 11 c ******************************************************************** 12 c Shema d'advection " pseudo amont " dans la verticale13 c pour appel dans la physique(sedimentation)12 c "pseudo upstream" Advection scheme along the vertical 13 c to be used in the physics (sedimentation) 14 14 c ******************************************************************** 15 c q rapport de melange (kg/kg)... 16 c masse : masse de la couche Dp/g 17 c w : masse d'atm ``transferee'' a chaque pas de temps (kg.m-2) 18 c pente_max = 2 conseillee 19 c 20 c 21 c -------------------------------------------------------------------- 15 22 16 IMPLICIT NONE 23 17 c … … 29 23 integer,intent(in) :: ngrid ! number of atmospheric columns 30 24 integer,intent(in) :: nlay ! number of atmospheric layers 31 real masse(ngrid,nlay),pente_max 32 REAL q(ngrid,nlay) 33 REAL w(ngrid,nlay) 34 REAL wq(ngrid,nlay+1) 25 real,intent(in) :: masse(ngrid,nlay) ! mass of atmospheric layer delta(P)/g 26 real,intent(in) :: pente_max ! maximum slope for the scheme (2 is recommended) 27 real,intent(inout) :: q(ngrid,nlay) ! tracer mixing ratio (kg/kg) 28 real,intent(inout) :: w(ngrid,nlay) ! mass of atmosphere "transfered" over the time step (kg.m-2) 29 real,intent(out) :: wq(ngrid,nlay+1) ! trancer increment due to advection (kg) 35 30 c 36 31 c Local … … 45 40 integer m 46 41 47 REAL SSUM,CVMGP,CVMGT48 integer ismax,ismin49 42 50 51 c On oriente tout dans le sens de la pression c'est a dire dans le 52 c sens de W 43 c Orientation follows pressure, i.e. follows W 53 44 54 45 do l=2,nlay … … 61 52 do l=2,nlay-1 62 53 do ij=1,ngrid 63 #ifdef CRAY64 dzq(ij,l)=0.5*65 , cvmgp(dzqw(ij,l)+dzqw(ij,l+1),0.,dzqw(ij,l)*dzqw(ij,l+1))66 #else67 54 if(dzqw(ij,l)*dzqw(ij,l+1).gt.0.) then 68 55 dzq(ij,l)=0.5*(dzqw(ij,l)+dzqw(ij,l+1)) … … 70 57 dzq(ij,l)=0. 71 58 endif 72 #endif 59 73 60 dzqmax=pente_max*min(adzqw(ij,l),adzqw(ij,l+1)) 74 61 dzq(ij,l)=sign(min(abs(dzq(ij,l)),dzqmax),dzq(ij,l)) … … 81 68 enddo 82 69 c --------------------------------------------------------------- 83 c .... c alcul des termes d'advection verticale.......70 c .... compute vertical advection terms ....... 84 71 c --------------------------------------------------------------- 85 72 86 c c alcul de - d( q * w )/ d(sigma) qu'on ajoute a dq pour calculerdq73 c compute - d( q * w )/ d(sigma) , later added to dq to compute dq 87 74 c 88 75 c No flux at the model top:
Note: See TracChangeset
for help on using the changeset viewer.