| 1 | module age_of_air_mod |
|---|
| 2 | |
|---|
| 3 | !**************************************************************** |
|---|
| 4 | ! |
|---|
| 5 | ! Age of air module |
|---|
| 6 | ! |
|---|
| 7 | ! Author: Maureen Cohen |
|---|
| 8 | ! ------- |
|---|
| 9 | ! |
|---|
| 10 | ! Version: 01/12/2023 |
|---|
| 11 | ! See "Age of air as a diagnostic for transport timescales in global models," |
|---|
| 12 | ! Geosci. Model Dev., 11, 3109–3130, 2018, https://doi.org/10.5194/gmd-11-3109-2018 |
|---|
| 13 | ! |
|---|
| 14 | !***************************************************************** |
|---|
| 15 | |
|---|
| 16 | implicit none |
|---|
| 17 | |
|---|
| 18 | logical,save :: ok_aoa ! flag to trigger runing with Age of Air scheme |
|---|
| 19 | logical,save :: reinit_aoa ! flag to reinitialize Age of air scheme |
|---|
| 20 | integer, save :: i_aoa ! index of age of air tracer in the tracer array |
|---|
| 21 | !$OMP THREADPRIVATE(ok_aoa,reinit_aoa,i_aoa) |
|---|
| 22 | real (kind=8), parameter :: forc_con = 1e-15 ! forcing constant in kg/kg/s |
|---|
| 23 | integer,save :: lev_aoa ! atmospheric level at which forcing is introduced |
|---|
| 24 | !$OMP THREADPRIVATE(lev_aoa) |
|---|
| 25 | real, allocatable, save:: init_age(:,:) ! saved initial value of age(:,:) |
|---|
| 26 | !$OMP THREADPRIVATE(init_age) |
|---|
| 27 | real, allocatable, save :: init_aoa(:,:) ! saved initial value of qx(:,:,i_aoa) |
|---|
| 28 | !$OMP THREADPRIVATE(init_aoa) |
|---|
| 29 | |
|---|
| 30 | !---------------------------------------------------------------------------- |
|---|
| 31 | contains |
|---|
| 32 | !---------------------------------------------------------------------------- |
|---|
| 33 | |
|---|
| 34 | subroutine age_of_air(aoa_tracer,klon,klev,itap,pdtphys,age_offset,trac_offset,aoa) |
|---|
| 35 | |
|---|
| 36 | implicit none |
|---|
| 37 | |
|---|
| 38 | ! Inputs |
|---|
| 39 | integer,intent(in) :: klon, klev ! number of gridpoints and levels |
|---|
| 40 | integer,intent(in) :: itap ! physics timestep counter |
|---|
| 41 | real,intent(in) :: pdtphys ! physics timestep length (s) |
|---|
| 42 | real,intent(in) :: age_offset(klon,klev) ! Init value of age, from restartphy or reinialitised to 0 |
|---|
| 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) |
|---|
| 45 | |
|---|
| 46 | ! In/out |
|---|
| 47 | real,intent(inout) :: aoa_tracer(klon,klev) ! age of air tracer mixing ratio (kg/kg) |
|---|
| 48 | |
|---|
| 49 | ! Outputs |
|---|
| 50 | real,intent(out) :: aoa(klon,klev) ! age of air value in seconds |
|---|
| 51 | |
|---|
| 52 | ! Local |
|---|
| 53 | real :: ts_forc ! forcing added in this timestep in kg/kg |
|---|
| 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 |
|---|
| 58 | |
|---|
| 59 | |
|---|
| 60 | !******************Add forcing constant to level lev_aoa********** |
|---|
| 61 | |
|---|
| 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 |
|---|
| 71 | |
|---|
| 72 | aoa_tracer(:,lev_aoa) = ts_forc |
|---|
| 73 | |
|---|
| 74 | aoa(:,:) = (run_time + cumul_time) - aoa_tracer(:,:)/forc_con ! calculate age of air in seconds |
|---|
| 75 | |
|---|
| 76 | end subroutine age_of_air |
|---|
| 77 | |
|---|
| 78 | !***************************************************************** |
|---|
| 79 | ! |
|---|
| 80 | ! Extracts index of age of air tracer in the tracer array |
|---|
| 81 | ! |
|---|
| 82 | !***************************************************************** |
|---|
| 83 | subroutine aoa_ini(ok_chem, tr_scheme) |
|---|
| 84 | use infotrac_phy, only: nqtot, tname |
|---|
| 85 | use chemparam_mod, only: M_tr, type_tr |
|---|
| 86 | implicit none |
|---|
| 87 | |
|---|
| 88 | integer,intent(in) :: tr_scheme |
|---|
| 89 | logical,intent(in) :: ok_chem |
|---|
| 90 | integer :: i |
|---|
| 91 | character(len=*),parameter :: rname="aoa_ini" |
|---|
| 92 | |
|---|
| 93 | if (.not. ok_chem .and. (tr_scheme .eq. 0)) then |
|---|
| 94 | allocate(M_tr(nqtot)) ! allocate M_tr only if this has not already been done in another routine |
|---|
| 95 | allocate(type_tr(nqtot)) ! same for type_tr |
|---|
| 96 | endif |
|---|
| 97 | |
|---|
| 98 | ! Initialise index of aoa to 0 |
|---|
| 99 | i_aoa = 0 |
|---|
| 100 | |
|---|
| 101 | do i=1, nqtot |
|---|
| 102 | |
|---|
| 103 | write(*,*) rname,' i',i |
|---|
| 104 | write(*,*) rname,' tname(i): ',tname(i) |
|---|
| 105 | |
|---|
| 106 | select case(tname(i)) |
|---|
| 107 | case('aoa') |
|---|
| 108 | i_aoa=i |
|---|
| 109 | write(*,*) rname,' aoa: ',i_aoa |
|---|
| 110 | M_tr(i_aoa) = 44.0 |
|---|
| 111 | write(*,*) rname,' M_tr(i_aoa): ', 44.0 |
|---|
| 112 | type_tr(i_aoa) = 10 |
|---|
| 113 | |
|---|
| 114 | end select |
|---|
| 115 | end do |
|---|
| 116 | |
|---|
| 117 | end subroutine aoa_ini |
|---|
| 118 | |
|---|
| 119 | end module age_of_air_mod |
|---|