Changeset 3622 for trunk/LMDZ.VENUS
- Timestamp:
- Feb 12, 2025, 5:16:29 PM (21 hours ago)
- Location:
- trunk/LMDZ.VENUS/libf/phyvenus
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.VENUS/libf/phyvenus/age_of_air_mod.F90
r3451 r3622 25 25 real, allocatable, save:: init_age(:,:) ! saved initial value of age(:,:) 26 26 !$OMP THREADPRIVATE(init_age) 27 real, allocatable, save :: init_aoa(:,:) ! saved initial value of qx(:,:,i_aoa) 28 !$OMP THREADPRIVATE(init_aoa) 27 29 28 30 !---------------------------------------------------------------------------- … … 30 32 !---------------------------------------------------------------------------- 31 33 32 subroutine age_of_air(aoa_tracer,klon,klev,itap,pdtphys,age_offset, aoa)34 subroutine age_of_air(aoa_tracer,klon,klev,itap,pdtphys,age_offset,trac_offset,aoa) 33 35 34 36 implicit none … … 40 42 real,intent(in) :: age_offset(klon,klev) ! Init value of age, from restartphy or reinialitised to 0 41 43 ! This value is the same every time step, but changes between restarts 44 real,intent(in) :: trac_offset(klon,klev) ! Same as above but for the tracer (kg/kg) rather than age (s) 42 45 43 46 ! In/out … … 49 52 ! Local 50 53 real :: ts_forc ! forcing added in this timestep in kg/kg 51 real :: run_time ! elapsed simulation time in seconds 54 real :: run_time ! elapsed simulation time in seconds (this run only) 55 real :: past_time(klon,1) ! time from past runs in seconds in source level only (if job chaining) 56 real :: cumul_time ! average of past_time array 57 integer :: i 52 58 53 59 54 60 !******************Add forcing constant to level lev_aoa********** 55 61 56 run_time = (itap+1)*pdtphys ! (number of physics timesteps * seconds per physics timestep) = seconds elapsed since start 57 58 ts_forc = run_time*forc_con ! elapsed time x forcing constant = tracer source this timestep 62 run_time = (itap+1)*pdtphys ! (number of physics timesteps * seconds per physics timestep) = seconds elapsed since start of this run 63 64 do i=1,klon 65 past_time = trac_offset(i,lev_aoa)/forc_con ! Convert initial tracer mixing ratio in source region into elapsed time from past runs 66 enddo 67 68 cumul_time = sum(past_time)/size(past_time) ! Average so we have a scalar to match run_time 69 70 ts_forc = (run_time + cumul_time)*forc_con ! total elapsed time x forcing constant = tracer source this timestep 59 71 60 72 aoa_tracer(:,lev_aoa) = ts_forc 61 73 62 aoa(:,:) = run_time - aoa_tracer(:,:)/forc_con + age_offset(:,:)! calculate age of air in seconds74 aoa(:,:) = (run_time + cumul_time) - aoa_tracer(:,:)/forc_con ! calculate age of air in seconds 63 75 64 76 end subroutine age_of_air -
trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F
r3615 r3622 69 69 USE chemparam_mod 70 70 USE age_of_air_mod, ONLY: ok_aoa, reinit_aoa, i_aoa, init_age 71 USE age_of_air_mod, ONLY: aoa_ini, age_of_air 71 USE age_of_air_mod, ONLY: aoa_ini, age_of_air, init_aoa 72 72 USE conc 73 73 USE param_v4_h … … 749 749 write(lunout,*) 'Initialising age of air subroutine' 750 750 allocate(init_age(klon,klev)) 751 allocate(init_aoa(klon,klev)) 751 752 call aoa_ini(ok_chem, tr_scheme) ! get index of age of air tracer in the tracer array 752 753 if (reinit_aoa) then 753 754 age(:,:) = 0. ! Set initial value of age of air to 0 everywhere 754 tr_seri(:,:,i_aoa) = 1e-30 ! Set initial value of tracer to tiny everywhere 755 init_aoa(:,:) = 1e-30 ! Set initial value of tracer to tiny everywhere 756 else 757 init_aoa(:,:) = qx(:,:,i_aoa) 755 758 endif ! reinitialisation loop 756 759 init_age(:,:) = age(:,:) ! save initial value of age, either read in from restartphy or set to 0 760 !init_aoa(:,:) = tr_seri(:,:,i_aoa) ! same for tracer 757 761 endif ! age of air initialisation 758 762 … … 814 818 ENDDO 815 819 ENDDO 820 821 ! Handle reinitialization of age of air tracer case 822 if (debut.and.ok_aoa.and.reinit_aoa) then 823 tr_seri(:,:,i_aoa)=init_aoa(:,:) 824 endif 816 825 C 817 826 DO i = 1, klon … … 1056 1065 if ( ok_aoa ) then 1057 1066 call age_of_air(tr_seri(:,:,i_aoa),klon,klev, 1058 I itap,pdtphys,init_age, 1067 I itap,pdtphys,init_age,init_aoa, 1059 1068 O age) 1060 1069 end if
Note: See TracChangeset
for help on using the changeset viewer.