Changeset 4151 for trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
- Timestamp:
- Mar 24, 2026, 5:50:52 PM (13 days ago)
- File:
-
- 1 edited
-
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F (modified) (4 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r4150 r4151 98 98 use nonoro_gwd_ran_mod, only: nonoro_gwd_ran 99 99 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 100 105 use check_fields_mod, only: check_physics_fields 101 106 use surfini_mod, only: surfini … … 584 589 REAL :: viscoh ! kinematic molecular viscosity for heat 585 590 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 586 595 c======================================================================= 587 596 pdq(:,:,:) = 0. … … 878 887 zplay(:,:) = pplay(:,:) 879 888 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) 880 901 881 902 !#ifndef MESOSCALE … … 2568 2589 ENDDO 2569 2590 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 2570 2596 2571 2597 ! Density
Note: See TracChangeset
for help on using the changeset viewer.
