Changeset 4043 for trunk/LMDZ.MARS/libf
- Timestamp:
- Feb 5, 2026, 2:47:26 PM (4 weeks ago)
- Location:
- trunk/LMDZ.MARS/libf/phymars/dyn1d
- Files:
-
- 3 edited
-
init_testphys1d_mod.F90 (modified) (3 diffs)
-
testphys1d.F90 (modified) (2 diffs)
-
writerestart1D_mod.F90 (modified) (2 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/dyn1d/init_testphys1d_mod.F90
r3996 r4043 313 313 ! ------------------------ 314 314 psurf = 610. ! Default value for psurf 315 pa = 20. ! Default value for transition pressure (for hybrid coord.) 316 preff = 610. ! Default value for reference surface pressure 315 317 write(*,*) 'Surface pressure (Pa)?' 316 318 if (.not. therestart1D) then … … 318 320 else 319 321 open(3,file = trim(start1Dname),status = "old",action = "read") 320 read(3,*) header, psurf 321 endif 322 write(*,*) " psurf = ",psurf 322 read(3,*) header, psurf, pa, preff 323 endif 324 write(*,*) " psurf = ", psurf 325 write(*,*) " pa = ", pa 326 write(*,*) " preff = ", preff 323 327 324 328 ! Compute pressures and altitudes of atmospheric levels … … 330 334 call getin("hybrid",hybrid) 331 335 write(*,*) " hybrid = ", hybrid 332 333 ! Reference pressures334 pa = 20. ! Transition pressure (for hybrid coord.)335 preff = 610. ! Reference surface pressure336 336 337 337 call disvert_noterre -
trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90
r3912 r4043 11 11 use dimradmars_mod, only: tauvis, totcloudfrac, albedo 12 12 use dust_param_mod, only: tauscaling 13 use comvert_mod, only: ap, bp, aps, bps 13 use comvert_mod, only: ap, bp, aps, bps, pa, preff 14 14 use physiq_mod, only: physiq 15 15 use turb_mod, only: q2 … … 261 261 262 262 ! Writing the "restart1D.txt" file for the next run 263 if (startfiles_1D) call writerestart1D('restart1D.txt',psurf, tsurf(1,:),nlayer,size(tsurf,2),temp,u,v,nq,noms,qsurf(1,:,:),q)263 if (startfiles_1D) call writerestart1D('restart1D.txt',psurf,pa,preff,tsurf(1,:),nlayer,size(tsurf,2),temp,u,v,nq,noms,qsurf(1,:,:),q) 264 264 265 265 write(*,*) "testphys1d: everything is cool!" -
trunk/LMDZ.MARS/libf/phymars/dyn1d/writerestart1D_mod.F90
r3074 r4043 5 5 contains 6 6 7 SUBROUTINE writerestart1D(filename,psurf, tsurf,nlayer,nslope,temp,u,v,nq,qnames,qsurf,q)7 SUBROUTINE writerestart1D(filename,psurf,pa,preff,tsurf,nlayer,nslope,temp,u,v,nq,qnames,qsurf,q) 8 8 9 9 implicit none 10 10 11 11 ! Arguments 12 character( len =*), intent(in) :: filename13 integer, intent(in) :: nlayer, nq, nslope14 real, intent(in) :: psurf15 real, dimension(nslope), intent(in) :: tsurf16 real, dimension(nlayer), intent(in) :: temp, u, v17 real, dimension(nlayer,nq), intent(in) :: q18 real, dimension(nq,nslope), intent(in) :: qsurf19 character( len =*), dimension(nq), intent(in) :: qnames12 character(*), intent(in) :: filename 13 integer, intent(in) :: nlayer, nq, nslope 14 real, intent(in) :: psurf, pa, preff 15 real, dimension(nslope), intent(in) :: tsurf 16 real, dimension(nlayer), intent(in) :: temp, u, v 17 real, dimension(nlayer,nq), intent(in) :: q 18 real, dimension(nq,nslope), intent(in) :: qsurf 19 character(*), dimension(nq), intent(in) :: qnames 20 20 21 21 ! Local variables … … 24 24 ! Write the data needed for a restart in "restart1D.txt" 25 25 open(1,file = filename,status = "replace",action = "write") 26 write(1,*) 'ps', psurf 26 write(1,*) 'ps', psurf, pa, preff 27 27 do i = 1,nq 28 28 write(1,*) qnames(i), (qsurf(i,j), j = 1,nslope), (q(il,i), il = 1,nlayer)
Note: See TracChangeset
for help on using the changeset viewer.
