Index: trunk/LMDZ.VENUS/libf/phyvenus/age_of_air_mod.F90
===================================================================
--- trunk/LMDZ.VENUS/libf/phyvenus/age_of_air_mod.F90	(revision 3615)
+++ trunk/LMDZ.VENUS/libf/phyvenus/age_of_air_mod.F90	(revision 3622)
@@ -25,4 +25,6 @@
 real, allocatable, save:: init_age(:,:) ! saved initial value of age(:,:)
 !$OMP THREADPRIVATE(init_age)
+real, allocatable, save :: init_aoa(:,:) ! saved initial value of qx(:,:,i_aoa)
+!$OMP THREADPRIVATE(init_aoa)
 
 !---------------------------------------------------------------------------- 
@@ -30,5 +32,5 @@
 !----------------------------------------------------------------------------
 
-subroutine age_of_air(aoa_tracer,klon,klev,itap,pdtphys,age_offset,aoa)
+subroutine age_of_air(aoa_tracer,klon,klev,itap,pdtphys,age_offset,trac_offset,aoa)
 
     implicit none
@@ -40,4 +42,5 @@
     real,intent(in) :: age_offset(klon,klev)   ! Init value of age, from restartphy or reinialitised to 0
                                                ! This value is the same every time step, but changes between restarts
+    real,intent(in) :: trac_offset(klon,klev)   ! Same as above but for the tracer (kg/kg) rather than age (s)
     
 ! In/out
@@ -49,16 +52,25 @@
 ! Local
     real :: ts_forc ! forcing added in this timestep in kg/kg
-    real :: run_time ! elapsed simulation time in seconds
+    real :: run_time ! elapsed simulation time in seconds (this run only)
+    real :: past_time(klon,1) ! time from past runs in seconds in source level only (if job chaining)
+    real :: cumul_time ! average of past_time array
+    integer :: i 
 
     
 !******************Add forcing constant to level lev_aoa**********
 
-    run_time = (itap+1)*pdtphys ! (number of physics timesteps * seconds per physics timestep) = seconds elapsed since start
-
-    ts_forc = run_time*forc_con ! elapsed time x forcing constant = tracer source this timestep
+    run_time = (itap+1)*pdtphys ! (number of physics timesteps * seconds per physics timestep) = seconds elapsed since start of this run
+    
+    do i=1,klon
+        past_time = trac_offset(i,lev_aoa)/forc_con ! Convert initial tracer mixing ratio in source region into elapsed time from past runs
+    enddo
+    
+    cumul_time = sum(past_time)/size(past_time) ! Average so we have a scalar to match run_time 
+    
+    ts_forc = (run_time + cumul_time)*forc_con ! total elapsed time x forcing constant = tracer source this timestep
     
     aoa_tracer(:,lev_aoa) = ts_forc
     
-    aoa(:,:) = run_time - aoa_tracer(:,:)/forc_con + age_offset(:,:) ! calculate age of air in seconds
+    aoa(:,:) = (run_time + cumul_time) - aoa_tracer(:,:)/forc_con ! calculate age of air in seconds
 
 end subroutine age_of_air
Index: trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F
===================================================================
--- trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F	(revision 3615)
+++ trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F	(revision 3622)
@@ -69,5 +69,5 @@
       USE chemparam_mod
       USE age_of_air_mod, ONLY: ok_aoa, reinit_aoa, i_aoa, init_age
-      USE age_of_air_mod, ONLY: aoa_ini, age_of_air
+      USE age_of_air_mod, ONLY: aoa_ini, age_of_air, init_aoa
       USE conc
       USE param_v4_h
@@ -749,10 +749,14 @@
           write(lunout,*) 'Initialising age of air subroutine'
           allocate(init_age(klon,klev))
+          allocate(init_aoa(klon,klev))
           call aoa_ini(ok_chem, tr_scheme) ! get index of age of air tracer in the tracer array 
           if (reinit_aoa) then
             age(:,:) = 0. ! Set initial value of age of air to 0 everywhere
-            tr_seri(:,:,i_aoa) = 1e-30 ! Set initial value of tracer to tiny everywhere
+            init_aoa(:,:) = 1e-30 ! Set initial value of tracer to tiny everywhere
+          else
+            init_aoa(:,:) = qx(:,:,i_aoa)
           endif ! reinitialisation loop
           init_age(:,:) = age(:,:) !  save initial value of age, either read in from restartphy or set to 0
+          !init_aoa(:,:) = tr_seri(:,:,i_aoa) ! same for tracer
       endif ! age of air initialisation
 
@@ -814,4 +818,9 @@
       ENDDO
       ENDDO
+
+      ! Handle reinitialization of age of air tracer case
+      if (debut.and.ok_aoa.and.reinit_aoa) then
+        tr_seri(:,:,i_aoa)=init_aoa(:,:)
+      endif
 C
       DO i = 1, klon
@@ -1056,5 +1065,5 @@
         if ( ok_aoa ) then
             call age_of_air(tr_seri(:,:,i_aoa),klon,klev,
-     I                   itap,pdtphys,init_age,
+     I                   itap,pdtphys,init_age,init_aoa,
      O                   age)
         end if
