[3451] | 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 | |
---|
| 28 | !---------------------------------------------------------------------------- |
---|
| 29 | contains |
---|
| 30 | !---------------------------------------------------------------------------- |
---|
| 31 | |
---|
| 32 | subroutine age_of_air(aoa_tracer,klon,klev,itap,pdtphys,age_offset,aoa) |
---|
| 33 | |
---|
| 34 | implicit none |
---|
| 35 | |
---|
| 36 | ! Inputs |
---|
| 37 | integer,intent(in) :: klon, klev ! number of gridpoints and levels |
---|
| 38 | integer,intent(in) :: itap ! physics timestep counter |
---|
| 39 | real,intent(in) :: pdtphys ! physics timestep length (s) |
---|
| 40 | real,intent(in) :: age_offset(klon,klev) ! Init value of age, from restartphy or reinialitised to 0 |
---|
| 41 | ! This value is the same every time step, but changes between restarts |
---|
| 42 | |
---|
| 43 | ! In/out |
---|
| 44 | real,intent(inout) :: aoa_tracer(klon,klev) ! age of air tracer mixing ratio (kg/kg) |
---|
| 45 | |
---|
| 46 | ! Outputs |
---|
| 47 | real,intent(out) :: aoa(klon,klev) ! age of air value in seconds |
---|
| 48 | |
---|
| 49 | ! Local |
---|
| 50 | real :: ts_forc ! forcing added in this timestep in kg/kg |
---|
| 51 | real :: run_time ! elapsed simulation time in seconds |
---|
| 52 | |
---|
| 53 | |
---|
| 54 | !******************Add forcing constant to level lev_aoa********** |
---|
| 55 | |
---|
| 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 |
---|
| 59 | |
---|
| 60 | aoa_tracer(:,lev_aoa) = ts_forc |
---|
| 61 | |
---|
| 62 | aoa(:,:) = run_time - aoa_tracer(:,:)/forc_con + age_offset(:,:) ! calculate age of air in seconds |
---|
| 63 | |
---|
| 64 | end subroutine age_of_air |
---|
| 65 | |
---|
| 66 | !***************************************************************** |
---|
| 67 | ! |
---|
| 68 | ! Extracts index of age of air tracer in the tracer array |
---|
| 69 | ! |
---|
| 70 | !***************************************************************** |
---|
| 71 | subroutine aoa_ini(ok_chem, tr_scheme) |
---|
| 72 | use infotrac_phy, only: nqtot, tname |
---|
| 73 | use chemparam_mod, only: M_tr, type_tr |
---|
| 74 | implicit none |
---|
| 75 | |
---|
| 76 | integer,intent(in) :: tr_scheme |
---|
| 77 | logical,intent(in) :: ok_chem |
---|
| 78 | integer :: i |
---|
| 79 | character(len=*),parameter :: rname="aoa_ini" |
---|
| 80 | |
---|
| 81 | if (.not. ok_chem .and. (tr_scheme .eq. 0)) then |
---|
| 82 | allocate(M_tr(nqtot)) ! allocate M_tr only if this has not already been done in another routine |
---|
| 83 | allocate(type_tr(nqtot)) ! same for type_tr |
---|
| 84 | endif |
---|
| 85 | |
---|
| 86 | ! Initialise index of aoa to 0 |
---|
| 87 | i_aoa = 0 |
---|
| 88 | |
---|
| 89 | do i=1, nqtot |
---|
| 90 | |
---|
| 91 | write(*,*) rname,' i',i |
---|
| 92 | write(*,*) rname,' tname(i): ',tname(i) |
---|
| 93 | |
---|
| 94 | select case(tname(i)) |
---|
| 95 | case('aoa') |
---|
| 96 | i_aoa=i |
---|
| 97 | write(*,*) rname,' aoa: ',i_aoa |
---|
| 98 | M_tr(i_aoa) = 44.0 |
---|
| 99 | write(*,*) rname,' M_tr(i_aoa): ', 44.0 |
---|
| 100 | type_tr(i_aoa) = 10 |
---|
| 101 | |
---|
| 102 | end select |
---|
| 103 | end do |
---|
| 104 | |
---|
| 105 | end subroutine aoa_ini |
---|
| 106 | |
---|
| 107 | end module age_of_air_mod |
---|