Changeset 3831
- Timestamp:
- Jul 4, 2025, 2:49:09 PM (18 hours ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/changelog.txt
r3824 r3831 4892 4892 == 02/07/2025 == JBC 4893 4893 Small corrections to the "display_netcdf.py" script in the "util" folder. 4894 4895 == 04/07/2025 == JBC 4896 Imposing the prescribed atmospheric water profile though 'dq' instead of 'q' directly (more correct). -
trunk/LMDZ.MARS/libf/phymars/dyn1d/init_testphys1d_mod.F90
r3741 r3831 790 790 if (prescribed_h2ovap) then 791 791 write(*,*) 'Atmospheric water profile is prescribed' 792 write(*,*) 'Unless it reaches saturation (maximal value)'793 792 open(10,file = 'profile_prescribed_h2o_vap',status = 'old',form = 'formatted',action = 'read',iostat = ierr) 794 793 if (ierr == 0) then … … 809 808 else 810 809 write(*,*) 'Atmospheric water profile is relaxed towards the profile in "profile_prescribed_h2o_vap"' 811 write(*,*) 'Unless it reaches saturation (maximal value)'812 810 endif 813 811 else -
trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90
r3741 r3831 6 6 use phyredem, only: physdem0, physdem1 7 7 use watersat_mod, only: watersat 8 use tracer_mod, only: igcm_h2o_vap, igcm_h2o_ice,noms8 use tracer_mod, only: igcm_h2o_vap, noms 9 9 use comcstfi_h, only: pi, g, rcp, cpp 10 10 use time_phylmdz_mod, only: daysec … … 188 188 189 189 ! Compute geopotential 190 ! ~~~~~~~~~~~~~~~~~~~~190 ! -------------------- 191 191 s = (aps/psurf + bps)**rcp 192 192 h = cpp*temp/(pks*s) … … 197 197 enddo 198 198 199 ! Prescription of atmospheric water if asked200 ! ------------------------------------------201 if (water) then202 call watersat(nlayer,temp,play,zqsat)203 if (prescribed_h2ovap) then ! If atmospheric water profile is prescribed204 if (h2ovap_relax_time < 0.) then ! Forcing205 q(1,:,igcm_h2o_vap) = min(zqsat(:),q_prescribed_h2o_vap(:))206 else ! Relaxation207 q(1,:,igcm_h2o_vap) = q_prescribed_h2o_vap(:) + (q(1,:,igcm_h2o_vap) - q_prescribed_h2o_vap(:))*dexp(-dttestphys/h2ovap_relax_time)208 endif209 endif210 endif211 212 199 ! Call physics 213 ! ------------ --------200 ! ------------ 214 201 call physiq(1,llm,nq,firstcall,lastcall,day,time,dttestphys,plev,play,phi,u,v,temp,q,w,du,dv,dtemp,dq,dpsurf) 215 202 ! ^----- outputs -----^ 216 203 217 204 ! Wind increment: specific for 1D 218 ! ------------------------------- -205 ! ------------------------------- 219 206 ! The physics compute the tendencies on u and v, 220 207 ! here we just add Coriolos effect … … 238 225 endif 239 226 240 ! Compute winds and temperaturefor next time step241 ! -------------------------------- ----------------227 ! Compute winds for next time step 228 ! -------------------------------- 242 229 u = u + dttestphys*du 243 230 v = v + dttestphys*dv 231 232 ! Compute temperature for next time step 233 ! -------------------------------------- 244 234 temp = temp + dttestphys*dtemp 245 235 … … 250 240 play = aps + psurf*bps 251 241 252 ! Increment tracers 242 ! Prescription of atmospheric water if asked 243 ! ------------------------------------------ 244 if (water .and. prescribed_h2ovap) then ! If atmospheric water profile is prescribed 245 if (h2ovap_relax_time < 0.) then ! Forcing 246 !call watersat(nlayer,temp,play,zqsat) 247 dq(1,:,igcm_h2o_vap) = (q_prescribed_h2o_vap(:) - q(1,:,igcm_h2o_vap))/dttestphys 248 else ! Relaxation 249 dq(1,:,igcm_h2o_vap) = (q_prescribed_h2o_vap(:) + (q(1,:,igcm_h2o_vap) - q_prescribed_h2o_vap(:))*dexp(-dttestphys/h2ovap_relax_time) - q(1,:,igcm_h2o_vap))/dttestphys 250 endif 251 endif 252 253 ! Compute tracers for next time step 254 ! ---------------------------------- 253 255 q = q + dttestphys*dq 254 256 enddo ! End of time stepping loop (idt=1,ndt)
Note: See TracChangeset
for help on using the changeset viewer.