Changeset 3912
- Timestamp:
- Sep 5, 2025, 2:46:19 PM (4 months ago)
- Location:
- trunk
- Files:
-
- 5 edited
-
LMDZ.COMMON/libf/evolution/pem.F90 (modified) (2 diffs)
-
LMDZ.MARS/changelog.txt (modified) (1 diff)
-
LMDZ.MARS/deftank/run.def.1d (modified) (1 diff)
-
LMDZ.MARS/libf/phymars/dyn1d/init_testphys1d_mod.F90 (modified) (3 diffs)
-
LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90 (modified) (4 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.COMMON/libf/evolution/pem.F90
r3907 r3912 276 276 277 277 ! Dummy variables to use the subroutine 'init_testphys1d' 278 logical :: therestart1D, therestartfi, prescribed_h2ovap278 logical :: therestart1D, therestartfi, ctrl_h2ovap, ctrl_h2oice 279 279 integer :: ndt, day0 280 real :: ptif, pks, day, gru, grv, h2ovap_relax_time280 real :: ptif, pks, day, gru, grv, relaxtime_h2ovap, relaxtime_h2oice 281 281 real, dimension(:), allocatable :: zqsat 282 282 real, dimension(:,:,:), allocatable :: dq, dqdyn 283 real, dimension(nlayer) :: play, w, q _prescribed_h2o_vap283 real, dimension(nlayer) :: play, w, qref_h2ovap, qref_h2oice 284 284 real, dimension(nlayer + 1) :: plev 285 285 #else … … 380 380 call init_testphys1d('start1D.txt','startfi.nc',therestart1D,therestartfi,ngrid,nlayer,610.,nq,q, & 381 381 time_0,pdyn(1),ucov,vcov,teta,ndt,ptif,pks,dtphys,zqsat,dq,dqdyn,day0,day,gru,grv,w, & 382 play,plev,latitude,longitude,cell_area,prescribed_h2ovap,h2ovap_relax_time,q_prescribed_h2o_vap) 382 play,plev,latitude,longitude,cell_area, & 383 ctrl_h2ovap,relaxtime_h2ovap,qref_h2ovap,ctrl_h2oice,relaxtime_h2oice,qref_h2oice) 383 384 nsplit_phys = 1 384 385 day_step = steps_per_sol -
trunk/LMDZ.MARS/changelog.txt
r3906 r3912 4968 4968 Fixing index bug in soilwater (index kept constant in a loop). 4969 4969 4970 == 05/09/2025 == JBC 4971 Prescription of the atmospheric water ice profile in 1D. The boolean 'ctrl_h2oice' activates the option, the file "profile_ref_h2ovap" defines the prescribed profile and the value 'relaxtime_h2oice', if positive, activates the relaxation and gives the constant of time. -
trunk/LMDZ.MARS/deftank/run.def.1d
r3741 r3912 49 49 # 0 == Northward, 90 == Eastward, 180 == Southward, 270 == Westward 50 50 slope_orientation=0. 51 # Prescription of atmospheric water profile defined in file "profile_prescribed_h2o_vap" (kg/kg)52 prescribed_h2ovap=.false.51 # Prescription of atmospheric water vapour profile defined in file "profile_ref_h2o_vap" (kg/kg) 52 ctrl_h2ovap=.false. 53 53 # Relaxation time (s). If <0 then no time relaxation i.e. forcing 54 h2ovap_relax_time=-1. 54 relaxtime_h2ovap=-1. 55 # Prescription of atmospheric water ice profile defined in file "profile_ref_h2o_ice" (kg/kg) 56 ctrl_h2oice=.false. 57 # Relaxation time (s). If <0 then no time relaxation i.e. forcing 58 relaxtime_h2oice=-1. 55 59 56 60 # hybrid vertical coordinate ? (.true. for hybrid and .false. for sigma levels) -
trunk/LMDZ.MARS/libf/phymars/dyn1d/init_testphys1d_mod.F90
r3855 r3912 7 7 SUBROUTINE init_testphys1d(start1Dname,startfiname,therestart1D,therestartfi,ngrid,nlayer,odpref, & 8 8 nq,q,time,psurf,u,v,temp,ndt,ptif,pks,dttestphys,zqsat,dq,dqdyn,day0,day,gru,grv,w, & 9 play,plev,latitude,longitude,cell_area,prescribed_h2ovap,h2ovap_relax_time,q_prescribed_h2o_vap) 9 play,plev,latitude,longitude,cell_area,ctrl_h2ovap,relaxtime_h2ovap,qref_h2ovap, & 10 ctrl_h2oice,relaxtime_h2oice,qref_h2oice) 10 11 11 12 use ioipsl_getincom, only: getin ! To use 'getin' … … 65 66 integer, intent(inout) :: nq 66 67 67 real, dimension(:,:,:), allocatable, intent(out) :: q ! tracer mixing ratio (e.g. kg/kg)68 real, intent(out) :: time ! time (0<time<1; time=0.5 at noon)69 real, intent(out) :: psurf ! Surface pressure70 real, dimension(nlayer), intent(out) :: u, v ! zonal, meridional wind71 real, dimension(nlayer), intent(out) :: temp ! temperature at the middle of the layers68 real, dimension(:,:,:), allocatable, intent(out) :: q ! tracer mixing ratio (e.g. kg/kg) 69 real, intent(out) :: time ! time (0<time<1; time=0.5 at noon) 70 real, intent(out) :: psurf ! Surface pressure 71 real, dimension(nlayer), intent(out) :: u, v ! zonal, meridional wind 72 real, dimension(nlayer), intent(out) :: temp ! temperature at the middle of the layers 72 73 integer, intent(out) :: ndt 73 74 real, intent(out) :: ptif, pks 74 real, intent(out) :: dttestphys ! testphys1d timestep75 real, intent(out) :: dttestphys ! testphys1d timestep 75 76 real, dimension(:), allocatable, intent(out) :: zqsat 76 real, dimension(:,:,:), allocatable, intent(out) :: dq, dqdyn ! Physical and dynamical tandencies77 integer, intent(out) :: day0 ! initial (sol ; =0 at Ls=0) and final date78 real, intent(out) :: day ! date during the run79 real, intent(out) :: gru, grv ! prescribed "geostrophic" background wind80 real, dimension(nlayer), intent(out) :: w ! "Dummy wind" in 1D81 real, dimension(nlayer), intent(out) :: play ! Pressure at the middle of the layers (Pa)82 real, dimension(nlayer + 1), intent(out) :: plev ! intermediate pressure levels (pa)77 real, dimension(:,:,:), allocatable, intent(out) :: dq, dqdyn ! Physical and dynamical tandencies 78 integer, intent(out) :: day0 ! initial (sol ; =0 at Ls=0) and final date 79 real, intent(out) :: day ! date during the run 80 real, intent(out) :: gru, grv ! prescribed "geostrophic" background wind 81 real, dimension(nlayer), intent(out) :: w ! "Dummy wind" in 1D 82 real, dimension(nlayer), intent(out) :: play ! Pressure at the middle of the layers (Pa) 83 real, dimension(nlayer + 1), intent(out) :: plev ! intermediate pressure levels (pa) 83 84 real, dimension(1), intent(out) :: latitude, longitude, cell_area 84 logical, intent(out) :: prescribed_h2ovap ! Prescribe atmospheric water profile85 real, intent(out) :: h2ovap_relax_time ! Relaxtion time towards the atmospheric water profile86 real, dimension(nlayer), intent(out) :: q _prescribed_h2o_vap ! User-defined atmospheric water profile85 logical, intent(out) :: ctrl_h2ovap, ctrl_h2oice ! Flags to prescribe atmospheric water vapour/ice profiles 86 real, intent(out) :: relaxtime_h2ovap, relaxtime_h2oice ! Relaxation time towards the atmospheric water vapour/ice profiles 87 real, dimension(nlayer), intent(out) :: qref_h2ovap, qref_h2oice ! User-defined atmospheric water vapour/ice profiles 87 88 88 89 !======================================================================= … … 781 782 if (.not. therestartfi .and. paleoclimate) perennial_co2ice = 0. 782 783 783 ! Prescription of atmospheric water profile784 ! ----------------------------------------- 785 prescribed_h2ovap = .false. ! Default: free atmospheric water profile 786 h2ovap_relax_time = -1. ! Default: no time relaxation787 q _prescribed_h2o_vap= 0. ! Default784 ! Prescription of atmospheric water vapour/ice profiles 785 ! ---------------------------------------------- 786 ctrl_h2ovap = .false. ; ctrl_h2oice = .false. ! Default: free atmospheric water vapour/ice profiles 787 relaxtime_h2ovap = -1. ; relaxtime_h2oice = -1. ! Default: no time relaxation 788 qref_h2ovap = 0. ; qref_h2oice = 0. ! Default 788 789 if (water) then 789 write(*,*)'Force atmospheric water profile (mixing ratio in kg/kg)?'790 call getin(' prescribed_h2ovap',prescribed_h2ovap)791 write(*,*) ' prescribed_h2ovap = ', prescribed_h2ovap792 if ( prescribed_h2ovap) then793 write(*,*) 'Atmospheric water profile is prescribed'794 open(10,file = 'profile_ prescribed_h2o_vap',status = 'old',form = 'formatted',action = 'read',iostat = ierr)790 write(*,*)'Force atmospheric water vapour profile (mixing ratio in kg/kg)?' 791 call getin('ctrl_h2ovap',ctrl_h2ovap) 792 write(*,*) 'ctrl_h2ovap = ', ctrl_h2ovap 793 if (ctrl_h2ovap) then 794 write(*,*) 'Atmospheric water vapour profile is prescribed' 795 open(10,file = 'profile_ref_h2o_vap',status = 'old',form = 'formatted',action = 'read',iostat = ierr) 795 796 if (ierr == 0) then 796 797 do ilayer = 1,nlayer 797 read(10,*,iostat=ierr) q _prescribed_h2o_vap(ilayer)798 if (ierr /= 0) error stop 'Not enough atmospheric layers defined in the file "profile_ prescribed_h2o_vap"!'798 read(10,*,iostat=ierr) qref_h2ovap(ilayer) 799 if (ierr /= 0) error stop 'Not enough atmospheric layers defined in the file "profile_ref_h2o_vap"!' 799 800 enddo 800 801 else 801 error stop 'File "profile_ prescribed_h2o_vap" was not found!'802 error stop 'File "profile_ref_h2o_vap" was not found!' 802 803 endif 803 804 close(10) 804 805 805 write(*,*) 'Relax atmospheric water profile with time constant (s)?'806 write(*,*) 'Relax atmospheric water vapour profile with time constant (s)?' 806 807 807 call getin(' h2ovap_relax_time',h2ovap_relax_time)808 write(*,*) ' h2ovap_relax_time = ', h2ovap_relax_time809 if ( h2ovap_relax_time< 0.) then810 write(*,*) 'Atmospheric water vapo r profile is not relaxed (forcing).'808 call getin('relaxtime_h2ovap',relaxtime_h2ovap) 809 write(*,*) 'relaxtime_h2ovap = ', relaxtime_h2ovap 810 if (relaxtime_h2ovap < 0.) then 811 write(*,*) 'Atmospheric water vapour profile is not relaxed (forcing).' 811 812 else 812 write(*,*) 'Atmospheric water profile is relaxed towards the profile in "profile_prescribed_h2o_vap"'813 write(*,*) 'Atmospheric water vapour profile is relaxed towards the profile in "profile_ref_h2o_vap"' 813 814 endif 814 815 else 815 write(*,*) 'Free atmospheric water vapor profile' 816 write(*,*) 'Total water is conserved in the column' 816 write(*,*) 'Free atmospheric water vapour profile' 817 endif 818 819 write(*,*)'Force atmospheric water ice profile (mixing ratio in kg/kg)?' 820 call getin('ctrl_h2oice',ctrl_h2oice) 821 write(*,*) 'ctrl_h2oice = ', ctrl_h2oice 822 if (ctrl_h2oice) then 823 write(*,*) 'Atmospheric water ice profile is prescribed' 824 open(10,file = 'profile_ref_h2o_ice',status = 'old',form = 'formatted',action = 'read',iostat = ierr) 825 if (ierr == 0) then 826 do ilayer = 1,nlayer 827 read(10,*,iostat=ierr) qref_h2oice(ilayer) 828 if (ierr /= 0) error stop 'Not enough atmospheric layers defined in the file "profile_ref_h2o_ice"!' 829 enddo 830 else 831 error stop 'File "profile_ref_h2o_ice" was not found!' 832 endif 833 close(10) 834 835 write(*,*) 'Relax atmospheric water ice profile with time constant (s)?' 836 837 call getin('relaxtime_h2oice',relaxtime_h2oice) 838 write(*,*) 'relaxtime_h2oice = ', relaxtime_h2oice 839 if (relaxtime_h2oice < 0.) then 840 write(*,*) 'Atmospheric water ice profile is not relaxed (forcing).' 841 else 842 write(*,*) 'Atmospheric water ice profile is relaxed towards the profile in "profile_ref_h2o_ice"' 843 endif 844 else 845 write(*,*) 'Free atmospheric water ice profile' 817 846 endif 818 847 endif -
trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90
r3838 r3912 6 6 use phyredem, only: physdem0, physdem1 7 7 use watersat_mod, only: watersat 8 use tracer_mod, only: igcm_h2o_vap, noms8 use tracer_mod, only: igcm_h2o_vap, igcm_h2o_ice, noms 9 9 use comcstfi_h, only: pi, g, rcp, cpp 10 10 use time_phylmdz_mod, only: daysec … … 83 83 logical :: startfiles_1D, therestart1D, therestartfi, there 84 84 85 ! Prescription of atmospheric water profile86 logical :: prescribed_h2ovap87 real :: h2ovap_relax_time88 real, dimension(nlayer) :: q _prescribed_h2o_vap85 ! Prescription of atmospheric water/ice profiles 86 logical :: ctrl_h2ovap, ctrl_h2oice 87 real :: relaxtime_h2ovap, relaxtime_h2oice 88 real, dimension(nlayer) :: qref_h2ovap, qref_h2oice 89 89 real, dimension(:), allocatable :: zqsat 90 90 !======================================================================= … … 142 142 call init_testphys1d('start1D.txt','startfi.nc',therestart1D,therestartfi,ngrid,nlayer,odpref,nq,q, & 143 143 time,psurf,u,v,temp,ndt,ptif,pks,dttestphys,zqsat,dq,dqdyn,day0,day,gru,grv,w, & 144 play,plev,latitude,longitude,cell_area,prescribed_h2ovap,h2ovap_relax_time,q_prescribed_h2o_vap) 144 play,plev,latitude,longitude,cell_area, & 145 ctrl_h2ovap,relaxtime_h2ovap,qref_h2ovap,ctrl_h2oice,relaxtime_h2oice,qref_h2oice) 145 146 146 147 ! Write a "startfi" file … … 208 209 ! Water tracer increment: specific for 1D 209 210 ! --------------------------------------- 210 ! The physics computes the tendency on q(1,:,igcm_h2o_vap), here we can mimic an effect of the "3D dynamics" 211 ! either by forcing a profile or by relaxing towards a prescribed profile 212 if (water .and. prescribed_h2ovap) then ! If atmospheric water profile is prescribed 213 if (h2ovap_relax_time < 0.) then ! Forcing 214 ! For some tests: unless it reaches saturation (maximal value) 215 !call watersat(nlayer,temp,play,zqsat) 216 !dq(1,:,igcm_h2o_vap) = (min(zqsat(:),q_prescribed_h2o_vap(:)) - q(1,:,igcm_h2o_vap))/dttestphys 217 dq(1,:,igcm_h2o_vap) = (q_prescribed_h2o_vap(:) - q(1,:,igcm_h2o_vap))/dttestphys 218 else ! Relaxation 219 dq(1,:,igcm_h2o_vap) = dq(1,:,igcm_h2o_vap) - (q(1,:,igcm_h2o_vap) - q_prescribed_h2o_vap(:))/h2ovap_relax_time 211 ! The physics computes the tendency on q for igcm_h2o_vap/igcm_h2o_ice, here we can mimic an effect 212 ! of the "3D dynamics" either by forcing a profile or by relaxing towards a prescribed profile 213 if (water) then 214 if (ctrl_h2ovap) then ! If atmospheric water vapour profile is prescribed 215 if (relaxtime_h2ovap < 0.) then ! Forcing 216 ! For some tests: unless it reaches saturation (maximal value) 217 !call watersat(nlayer,temp,play,zqsat) 218 !dq(1,:,igcm_h2o_vap) = (min(zqsat(:),qref_h2ovap(:)) - q(1,:,igcm_h2o_vap))/dttestphys 219 dq(1,:,igcm_h2o_vap) = (qref_h2ovap(:) - q(1,:,igcm_h2o_vap))/dttestphys 220 else ! Relaxation 221 dq(1,:,igcm_h2o_vap) = dq(1,:,igcm_h2o_vap) - (q(1,:,igcm_h2o_vap) - qref_h2ovap(:))/relaxtime_h2ovap 222 endif 223 endif 224 if (ctrl_h2oice) then ! If atmospheric water ice profile is prescribed 225 if (relaxtime_h2oice < 0.) then ! Forcing 226 dq(1,:,igcm_h2o_ice) = (qref_h2oice(:) - q(1,:,igcm_h2o_ice))/dttestphys 227 else ! Relaxation 228 dq(1,:,igcm_h2o_ice) = dq(1,:,igcm_h2o_ice) - (q(1,:,igcm_h2o_ice) - qref_h2oice(:))/relaxtime_h2oice 229 endif 220 230 endif 221 231 endif
Note: See TracChangeset
for help on using the changeset viewer.
