Changeset 3831


Ignore:
Timestamp:
Jul 4, 2025, 2:49:09 PM (18 hours ago)
Author:
jbclement
Message:

Mars PCM:
Imposing the prescribed atmospheric water profile though 'dq' instead of 'q' directly (more correct).
JBC

Location:
trunk/LMDZ.MARS
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/changelog.txt

    r3824 r3831  
    48924892== 02/07/2025 == JBC
    48934893Small corrections to the "display_netcdf.py" script in the "util" folder.
     4894
     4895== 04/07/2025 == JBC
     4896Imposing 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  
    790790    if (prescribed_h2ovap) then
    791791        write(*,*) 'Atmospheric water profile is prescribed'
    792         write(*,*) 'Unless it reaches saturation (maximal value)'
    793792        open(10,file = 'profile_prescribed_h2o_vap',status = 'old',form = 'formatted',action = 'read',iostat = ierr)
    794793        if (ierr == 0) then
     
    809808        else
    810809            write(*,*) 'Atmospheric water profile is relaxed towards the profile in "profile_prescribed_h2o_vap"'
    811             write(*,*) 'Unless it reaches saturation (maximal value)'
    812810        endif
    813811    else
  • trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F90

    r3741 r3831  
    66use phyredem,            only: physdem0, physdem1
    77use watersat_mod,        only: watersat
    8 use tracer_mod,          only: igcm_h2o_vap, igcm_h2o_ice, noms
     8use tracer_mod,          only: igcm_h2o_vap, noms
    99use comcstfi_h,          only: pi, g, rcp, cpp
    1010use time_phylmdz_mod,    only: daysec
     
    188188
    189189    ! Compute geopotential
    190     ! ~~~~~~~~~~~~~~~~~~~~
     190    ! --------------------
    191191    s = (aps/psurf + bps)**rcp
    192192    h = cpp*temp/(pks*s)
     
    197197    enddo
    198198
    199     ! Prescription of atmospheric water if asked
    200     ! ------------------------------------------
    201     if (water) then
    202         call watersat(nlayer,temp,play,zqsat)
    203         if (prescribed_h2ovap) then ! If atmospheric water profile is prescribed
    204             if (h2ovap_relax_time < 0.) then ! Forcing
    205                 q(1,:,igcm_h2o_vap) = min(zqsat(:),q_prescribed_h2o_vap(:))
    206             else ! Relaxation
    207                 q(1,:,igcm_h2o_vap) = q_prescribed_h2o_vap(:) + (q(1,:,igcm_h2o_vap) - q_prescribed_h2o_vap(:))*dexp(-dttestphys/h2ovap_relax_time)
    208             endif
    209         endif
    210      endif
    211 
    212199    ! Call physics
    213     ! --------------------
     200    ! ------------
    214201    call physiq(1,llm,nq,firstcall,lastcall,day,time,dttestphys,plev,play,phi,u,v,temp,q,w,du,dv,dtemp,dq,dpsurf)
    215202    !                                                                                      ^----- outputs -----^
    216203
    217204    ! Wind increment: specific for 1D
    218     ! --------------------------------
     205    ! -------------------------------
    219206    ! The physics compute the tendencies on u and v,
    220207    ! here we just add Coriolos effect
     
    238225    endif
    239226
    240     ! Compute winds and temperature for next time step
    241     ! ------------------------------------------------
     227    ! Compute winds for next time step
     228    ! --------------------------------
    242229    u = u + dttestphys*du
    243230    v = v + dttestphys*dv
     231
     232    ! Compute temperature for next time step
     233    ! --------------------------------------
    244234    temp = temp + dttestphys*dtemp
    245235
     
    250240    play = aps + psurf*bps
    251241
    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    ! ----------------------------------
    253255    q = q + dttestphys*dq
    254256enddo ! End of time stepping loop (idt=1,ndt)
Note: See TracChangeset for help on using the changeset viewer.