Changeset 3902
- Timestamp:
- Aug 21, 2025, 2:57:55 PM (4 months ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 15 edited
-
changelog.txt (modified) (1 diff)
-
libf/dynphy_lonlat/phymars/lect_start_archive.F (modified) (1 diff)
-
libf/dynphy_lonlat/phymars/newstart.F (modified) (1 diff)
-
libf/phymars/co2cloud.F90 (modified) (1 diff)
-
libf/phymars/geticecover.F90 (modified) (4 diffs)
-
libf/phymars/improvedco2clouds_mod.F90 (modified) (1 diff)
-
libf/phymars/interp_line.F (modified) (3 diffs)
-
libf/phymars/orbite.F (modified) (2 diffs)
-
libf/phymars/phyetat0_mod.F90 (modified) (1 diff)
-
libf/phymars/physiq_mod.F (modified) (2 diffs)
-
libf/phymars/simpleclouds.F (modified) (7 diffs)
-
libf/phymars/soil_settings.F (modified) (1 diff)
-
libf/phymars/tabfi.F (modified) (3 diffs)
-
libf/phymars/tcondco2.F90 (modified) (2 diffs)
-
libf/phymars/watercloud_mod.F (modified) (1 diff)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/changelog.txt
r3901 r3902 4951 4951 - clean module vlz_fi (remove obsolete #ifdef CRAY cpp alternatives) 4952 4952 - remove function cvmgp since it was only called under #ifdef CRAY in vlz_fi 4953 4954 == 21/08/2025 == EM 4955 More code tidying: turn geticecover, interp_line, orbite, simpleclouds, tabfi 4956 and tcondco2 into modules. -
trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/lect_start_archive.F
r3025 r3902 31 31 use comslope_mod, ONLY: nslope,def_slope,def_slope_mean, 32 32 & subslope_dist,end_comslope_h,ini_comslope_h 33 use interp_line_mod, only: interp_line 33 34 use netcdf 34 35 use surfdat_h, ONLY: watercaptag -
trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/newstart.F
r3886 r3902 38 38 use geometry_mod, only: longitude,latitude,cell_area 39 39 use lect_start_archive_mod, only: lect_start_archive 40 use tabfi_mod, only: tabfi 40 41 use phyetat0_mod, only: phyetat0 41 42 use phyredem, only: physdem0, physdem1 -
trunk/LMDZ.MARS/libf/phymars/co2cloud.F90
r3739 r3902 106 106 use datafile_mod, only: datadir 107 107 use density_co2_ice_mod, only: density_co2_ice 108 108 use tcondco2_mod, only: tcondco2 109 109 use improvedCO2clouds_mod, only: improvedCO2clouds 110 110 use microphys_h, only: nbinco2_cld, rad_cldco2, mco2 -
trunk/LMDZ.MARS/libf/phymars/geticecover.F90
r332 r3902 1 MODULE geticecover_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 SUBROUTINE geticecover( ngrid, ls, lontab, lattab, icecover ) 2 8 !!**********************************************!! … … 11 17 REAL, DIMENSION(ngrid), INTENT(OUT) :: icecover 12 18 !! LOCAL 13 REAL :: isitice14 19 INTEGER :: ig 15 20 icecover(:) = 0. … … 25 30 IMPLICIT NONE 26 31 REAL, INTENT(IN) :: ls,lon,lat 27 REAL :: nplatcrocus, splatcrocus28 32 isitice = 0. 29 33 !! for speedup purposes: … … 308 312 END FUNCTION splatcrocus 309 313 314 END MODULE geticecover_mod 315 310 316 !PROGRAM main 311 317 ! -
trunk/LMDZ.MARS/libf/phymars/improvedco2clouds_mod.F90
r3739 r3902 73 73 74 74 use conc_mod, only: mmean 75 use tcondco2_mod, only: tcondco2 75 76 use time_phylmdz_mod, only: daysec 76 77 use nucleaco2_mod, only: nucleaco2 -
trunk/LMDZ.MARS/libf/phymars/interp_line.F
r38 r3902 1 module interp_line_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine interp_line(x1,y1,len1,x2,y2,len2) 2 8 implicit none … … 14 20 ! --------- 15 21 ! inputs: 16 real x1(len1) ! ordered list of abscissas17 real y1(len1) ! values at x1(:)18 integer len1 ! length of x1(:) and y1(:)19 real x2(len2) !ordered list of abscissas at which interpolation is done20 integer len2 ! length of x2(:) and y2(:)22 real,intent(in) :: x1(len1) ! ordered list of abscissas 23 real,intent(in) :: y1(len1) ! values at x1(:) 24 integer,intent(in) :: len1 ! length of x1(:) and y1(:) 25 real,intent(in) :: x2(len2) !ordered list of abscissas at which interpolation is done 26 integer,intent(in) :: len2 ! length of x2(:) and y2(:) 21 27 ! outputs: 22 real y2(len2) ! interpolated values28 real,intent(out) :: y2(len2) ! interpolated values 23 29 !----------------------------------------------------------------------- 24 30 … … 53 59 enddo 54 60 55 end 61 end subroutine interp_line 62 63 end module interp_line_mod -
trunk/LMDZ.MARS/libf/phymars/orbite.F
r3040 r3902 1 MODULE orbite_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 SUBROUTINE orbite(pls,pdist_sol,pdecli) 2 8 USE planete_h, ONLY: e_elips, p_elips, obliquit, lsperi … … 37 43 pdecli = asin(sin(pls)*sin(obliquit*pi/180.)) 38 44 39 END 45 END SUBROUTINE orbite 46 47 END MODULE orbite_mod -
trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90
r3798 r3902 35 35 use geometry_mod, only: latitude 36 36 use soil_settings_mod, only: soil_settings 37 use tabfi_mod, only: tabfi 37 38 use callkeys_mod, only: startphy_file, rdstorm, hdo 38 39 -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r3901 r3902 64 64 use solang_mod, only: solang 65 65 use mucorr_mod, only: mucorr 66 use geticecover_mod, only: geticecover 66 67 use nirdata_mod, only: NIR_leedat 67 68 use nirco2abs_mod, only: nirco2abs … … 87 88 & obliquit 88 89 use planete_h, only: iniorbit 90 use orbite_mod, only: orbite 89 91 USE comcstfi_h, only: r, cpp, mugaz, g, rcp, pi, rad 90 92 USE calldrag_noro_mod, ONLY: calldrag_noro -
trunk/LMDZ.MARS/libf/phymars/simpleclouds.F
r3726 r3902 1 module simpleclouds_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine simpleclouds(ngrid,nlay,ptimestep, 2 8 & pplay,pzlay,pt,pdt, … … 22 28 c be found in the routine called "improvedclouds.F". 23 29 24 c Modif de zq si saturation dans l'atmosphere25 c si zq(ig,l)> zqsat(ig,l) -> zq(ig,l)=zqsat(ig,l)26 c Le test est effectue de bas en haut. L'eau condensee27 c (si saturation) est remise dans la couche en dessous.28 c L'eau condensee dans la couche du bas est deposee a la surface29 30 30 c Authors: Franck Montmessin (water ice scheme) 31 31 c Francois Forget (changed nuclei density & outputs) … … 39 39 c --------- 40 40 c Inputs: 41 INTEGER ngrid,nlay 42 integer nq ! nombre de traceurs 43 REAL ptimestep ! pas de temps physique (s) 44 REAL pplay(ngrid,nlay) ! pression au milieu des couches (Pa) 45 REAL pzlay(ngrid,nlay) ! altitude at the middle of the layers 46 REAL pt(ngrid,nlay) ! temperature at the middle of the 47 ! layers (K) 48 REAL pdt(ngrid,nlay) ! tendance temperature des autres 49 ! param. 50 real pq(ngrid,nlay,nq) ! traceur (kg/kg) 51 real pdq(ngrid,nlay,nq) ! tendance avant condensation 52 ! (kg/kg.s-1) 53 REAL tau(ngrid,naerkind) ! Column dust optical depth at each point 41 integer,intent(in) :: ngrid ! number of atmospheric columns 42 integer,intent(in) :: nlay ! number of atmospheric layers 43 integer,intent(in) :: nq ! number of tracers 44 real,intent(in) :: ptimestep ! physics time step (s) 45 real,intent(in) :: pplay(ngrid,nlay) ! pressure at mid-layer (Pa) 46 real,intent(in) :: pzlay(ngrid,nlay) ! altitude of the layers (m) 47 real,intent(in) :: pt(ngrid,nlay) ! input temperature at mid-layer (K) 48 real,intent(in) :: pdt(ngrid,nlay) ! tendency on temperature from previous paramatrizations (K/s) 49 real,intent(in) :: pq(ngrid,nlay,nq) ! input tracer mixing ratio (kg/kg) 50 real,intent(in) :: pdq(ngrid,nlay,nq) ! tendency on tracers from previous parametrizations (kg/kg.s-1) 51 real,intent(in) :: tau(ngrid,naerkind) ! Column dust optical depth in each column 54 52 55 53 c Output: 56 REAL rice(ngrid,nlay) ! Ice mass mean radius (m)57 ! (r_c in montmessin_2004)58 real pdqcloud(ngrid,nlay,nq) ! tendance de lacondensation59 ! H2O(kg/kg.s-1)60 REAL pdtcloud(ngrid,nlay) ! tendancetemperature due61 ! a la chaleur latente54 real,intent(out) :: rice(ngrid,nlay) ! Water ice mass mean radius (m) 55 ! (r_c in montmessin_2004) 56 real,intent(out) :: pdqcloud(ngrid,nlay,nq) ! tendencies due to H2O condensation 57 ! and sublimation (kg/kg.s-1) 58 real,intent(out) :: pdtcloud(ngrid,nlay) ! tendency on temperature due 59 ! to latent heat (K/s) 62 60 63 61 c------------------------------------------------------------------ … … 88 86 c ----------------- 89 87 90 c On "update" la valeur de q(nq) (water vapor) et temperature. 91 c On effectue qqes calculs preliminaires sur les couches : 88 c Update values of water vapor and ice, and temperature. 92 89 93 90 do l=1,nlay … … 126 123 c ---------------------------------------------- 127 124 c 128 c Rapport de melange a saturation dans la couche l : -------125 c Compute saturation mixing ratio in each layer 129 126 c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 130 127 131 128 call watersat(ngrid*nlay,zt,pplay,zqsat) 132 129 133 c taux de condensation (kg/kg/s-1) dans les differentes couches130 c compute condensation rates (kg/kg/s-1) in each layer 134 131 c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 135 132 … … 154 151 155 152 156 c Tendance finale153 c Final tendency 157 154 c ~~~~~~~~~~~~~~~ 158 155 do l=1, nlay … … 231 228 c endif !hdo 232 229 c------------------------------------------------------------------ 233 end 230 end subroutine simpleclouds 231 232 end module simpleclouds_mod -
trunk/LMDZ.MARS/libf/phymars/soil_settings.F
r3727 r3902 16 16 & inquire_field, inquire_dimension_length 17 17 USE comslope_mod, ONLY: nslope 18 use interp_line_mod, only: interp_line 18 19 implicit none 19 20 -
trunk/LMDZ.MARS/libf/phymars/tabfi.F
r3037 r3902 1 MODULE tabfi_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 c======================================================================= 2 8 SUBROUTINE tabfi(nid,Lmodif,tab0,day_ini,lmax,p_rad, … … 550 556 551 557 c----------------------------------------------------------------------- 552 end558 END SUBROUTINE tabfi 553 559 554 560 … … 584 590 585 591 586 592 END MODULE tabfi_mod 593 -
trunk/LMDZ.MARS/libf/phymars/tcondco2.F90
r2456 r3902 1 MODULE tcondco2_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 SUBROUTINE tcondco2(ngrid,nlay,p,q,tcond) 2 USE comcstfi_h 8 3 9 use conc_mod, only: mmean 4 10 … … 27 33 enddo 28 34 29 end 35 END SUBROUTINE tcondco2 36 37 END MODULE tcondco2_mod -
trunk/LMDZ.MARS/libf/phymars/watercloud_mod.F
r3871 r3902 18 18 USE updaterad, ONLY: updaterdust, updaterice_micro, 19 19 & updaterice_typ 20 USE simpleclouds_mod, ONLY: simpleclouds 20 21 USE improvedclouds_mod, ONLY: improvedclouds 21 22 USE watersat_mod, ONLY: watersat
Note: See TracChangeset
for help on using the changeset viewer.
