Changeset 4151 for trunk


Ignore:
Timestamp:
Mar 24, 2026, 5:50:52 PM (13 days ago)
Author:
yluo
Message:

Mars PCM:
Add missing updates to physiq_mod.F.
The previous commit only reverted the changes introduced in r4142 and omitted these updates.
YCL

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r4150 r4151  
    9898      use nonoro_gwd_ran_mod, only: nonoro_gwd_ran
    9999      use nonoro_gwd_mix_mod, only: nonoro_gwd_mix, calljliu_gwimix
     100      use tracer_mass_fixer_dyn_mod, only: call_mass_fixer_dyn,
     101     &                                   mass_predyn,
     102     &                                   tracer_mass_fixer_dyn,
     103     &                                   compute_tracer_mass_global,
     104     &                                   set_mass_predyn_from_startfi
    100105      use check_fields_mod, only: check_physics_fields
    101106      use surfini_mod, only: surfini
     
    584589      REAL :: viscoh                              ! kinematic molecular viscosity for heat
    585590
     591! Variables for correcting tracer mass non-conservation in dynamics
     592      real pdq_corrdyn(ngrid, nlayer, nq) ! tracer mass mixing ratio tendencies in order to correct
     593                                          ! mass non-conservations in dynamics
     594
    586595c=======================================================================
    587596      pdq(:,:,:) = 0.
     
    878887      zplay(:,:) = pplay(:,:)
    879888      ps(:) = pplev(:,1)
     889
     890! Correct mass non-conservation in dynamics
     891! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
     892      if (call_mass_fixer_dyn) then
     893         if (firstcall) then
     894            call set_mass_predyn_from_startfi
     895         endif ! if (firstcall)
     896
     897         call tracer_mass_fixer_dyn(ngrid, nlayer, nq, ptimestep,
     898     &                              pq, zplev, firstcall, pdq_corrdyn)
     899         pdq = pdq + pdq_corrdyn
     900      endif ! if (call_mass_fixer_dyn)
    880901
    881902!#ifndef MESOSCALE
     
    25682589        ENDDO
    25692590      ENDDO
     2591
     2592      if (call_mass_fixer_dyn) then     
     2593          call compute_tracer_mass_global(ngrid, nlayer, nq, zq,
     2594     &                                    zplev, mass_predyn)
     2595      end if
    25702596
    25712597      ! Density
Note: See TracChangeset for help on using the changeset viewer.