Changeset 3912 for trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90
- Timestamp:
- Sep 5, 2025, 2:46:19 PM (4 months ago)
- File:
-
- 1 edited
-
trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90 (modified) (4 diffs)
Legend:
- Unmodified
- Added
- Removed
-
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.
