| 1 | MODULE tracers |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! tracers |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! Tracer species management for tracking. |
|---|
| 8 | ! |
|---|
| 9 | ! AUTHORS & DATE |
|---|
| 10 | ! JB Clement, 12/2025 |
|---|
| 11 | ! |
|---|
| 12 | ! NOTES |
|---|
| 13 | ! |
|---|
| 14 | !----------------------------------------------------------------------- |
|---|
| 15 | |
|---|
| 16 | ! DEPENDENCIES |
|---|
| 17 | ! ------------ |
|---|
| 18 | use numerics, only: dp, di, k4 |
|---|
| 19 | |
|---|
| 20 | ! DECLARATION |
|---|
| 21 | ! ----------- |
|---|
| 22 | implicit none |
|---|
| 23 | |
|---|
| 24 | ! PARAMETERS |
|---|
| 25 | ! ---------- |
|---|
| 26 | character(11), parameter :: traceurdef_name = 'traceur.def' |
|---|
| 27 | integer(di), protected :: nq ! Number of tracers |
|---|
| 28 | character(30), dimension(:), allocatable, protected :: qnames ! Names of tracers |
|---|
| 29 | real(dp), dimension(:,:,:), allocatable, protected :: q_PCM ! Tracers in the "start.nc" at the beginning |
|---|
| 30 | integer(di), protected :: iPCM_qh2o ! Index for H2O vapor tracer from PCM |
|---|
| 31 | real(dp), dimension(:), allocatable, protected :: mmol ! Molar masses of tracers [g/mol] |
|---|
| 32 | |
|---|
| 33 | contains |
|---|
| 34 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 35 | |
|---|
| 36 | !======================================================================= |
|---|
| 37 | SUBROUTINE ini_tracers() |
|---|
| 38 | !----------------------------------------------------------------------- |
|---|
| 39 | ! NAME |
|---|
| 40 | ! ini_tracers |
|---|
| 41 | ! |
|---|
| 42 | ! DESCRIPTION |
|---|
| 43 | ! Allocate tracer molar mass array. |
|---|
| 44 | ! |
|---|
| 45 | ! AUTHORS & DATE |
|---|
| 46 | ! JB Clement, 12/2025 |
|---|
| 47 | ! |
|---|
| 48 | ! NOTES |
|---|
| 49 | ! |
|---|
| 50 | !----------------------------------------------------------------------- |
|---|
| 51 | |
|---|
| 52 | ! DEPENDENCIES |
|---|
| 53 | ! ------------ |
|---|
| 54 | use stoppage, only: stop_clean |
|---|
| 55 | use geometry, only: ngrid, nlayer |
|---|
| 56 | use display, only: print_msg |
|---|
| 57 | use utility, only: int2str |
|---|
| 58 | |
|---|
| 59 | ! DECLARATION |
|---|
| 60 | ! ----------- |
|---|
| 61 | implicit none |
|---|
| 62 | |
|---|
| 63 | ! LOCAL VARIABLES |
|---|
| 64 | ! --------------- |
|---|
| 65 | logical(k4) :: here |
|---|
| 66 | integer(di) :: i, ierr, funit |
|---|
| 67 | |
|---|
| 68 | ! CODE |
|---|
| 69 | ! ---- |
|---|
| 70 | inquire(file = traceurdef_name,exist = here) |
|---|
| 71 | if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//traceurdef_name//'"!',1) |
|---|
| 72 | call print_msg('> Reading "'//traceurdef_name//'"') |
|---|
| 73 | open(newunit = funit,file = traceurdef_name,status = 'old',form = 'formatted',action = 'read',iostat = ierr) |
|---|
| 74 | if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "'//traceurdef_name//'"!',ierr) |
|---|
| 75 | read(funit,*) nq |
|---|
| 76 | call print_msg('nq = '//int2str(nq)) |
|---|
| 77 | |
|---|
| 78 | ! Allocation |
|---|
| 79 | if (.not. allocated(mmol)) allocate(mmol(nq)) |
|---|
| 80 | if (.not. allocated(qnames)) allocate(qnames(nq)) |
|---|
| 81 | if (.not. allocated(q_PCM)) allocate(q_PCM(ngrid,nlayer,nq)) |
|---|
| 82 | |
|---|
| 83 | ! Initialize the names of tracers |
|---|
| 84 | do i = 1,nq |
|---|
| 85 | read(funit,*,iostat = ierr) qnames(i) |
|---|
| 86 | if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error reading the names of tracers in "'//traceurdef_name//'"!',1) |
|---|
| 87 | end do |
|---|
| 88 | close(funit) |
|---|
| 89 | |
|---|
| 90 | ! Initialize the indices of tracers |
|---|
| 91 | iPCM_qh2o = -1 |
|---|
| 92 | do i = 1,nq |
|---|
| 93 | if (qnames(i) == "h2o_vap") then |
|---|
| 94 | iPCM_qh2o = i |
|---|
| 95 | mmol(i) = 18._dp |
|---|
| 96 | end if |
|---|
| 97 | end do |
|---|
| 98 | |
|---|
| 99 | ! Checking if everything has been found |
|---|
| 100 | if (iPCM_qh2o < 0) call stop_clean(__FILE__,__LINE__,'H2O vapour index not found!',1) |
|---|
| 101 | |
|---|
| 102 | END SUBROUTINE ini_tracers |
|---|
| 103 | !======================================================================= |
|---|
| 104 | |
|---|
| 105 | !======================================================================= |
|---|
| 106 | SUBROUTINE end_tracers() |
|---|
| 107 | !----------------------------------------------------------------------- |
|---|
| 108 | ! NAME |
|---|
| 109 | ! end_tracers |
|---|
| 110 | ! |
|---|
| 111 | ! DESCRIPTION |
|---|
| 112 | ! Deallocate tracer molar mass array. |
|---|
| 113 | ! |
|---|
| 114 | ! AUTHORS & DATE |
|---|
| 115 | ! JB Clement, 12/2025 |
|---|
| 116 | ! |
|---|
| 117 | ! NOTES |
|---|
| 118 | ! |
|---|
| 119 | !----------------------------------------------------------------------- |
|---|
| 120 | |
|---|
| 121 | ! DECLARATION |
|---|
| 122 | ! ----------- |
|---|
| 123 | implicit none |
|---|
| 124 | |
|---|
| 125 | ! CODE |
|---|
| 126 | ! ---- |
|---|
| 127 | if (allocated(mmol)) deallocate(mmol) |
|---|
| 128 | if (allocated(qnames)) deallocate(qnames) |
|---|
| 129 | if (allocated(q_PCM)) deallocate(q_PCM) |
|---|
| 130 | |
|---|
| 131 | END SUBROUTINE end_tracers |
|---|
| 132 | !======================================================================= |
|---|
| 133 | |
|---|
| 134 | !======================================================================= |
|---|
| 135 | SUBROUTINE set_q_PCM(q_PCM_in,iq) |
|---|
| 136 | !----------------------------------------------------------------------- |
|---|
| 137 | ! NAME |
|---|
| 138 | ! set_q_PCM |
|---|
| 139 | ! |
|---|
| 140 | ! DESCRIPTION |
|---|
| 141 | ! Setter for 'q_PCM'. |
|---|
| 142 | ! |
|---|
| 143 | ! AUTHORS & DATE |
|---|
| 144 | ! JB Clement, 12/2025 |
|---|
| 145 | ! |
|---|
| 146 | ! NOTES |
|---|
| 147 | ! |
|---|
| 148 | !----------------------------------------------------------------------- |
|---|
| 149 | |
|---|
| 150 | ! DEPENDENCIES |
|---|
| 151 | ! ------------ |
|---|
| 152 | use geometry, only: dyngrd2vect, nlon, nlat, ngrid |
|---|
| 153 | |
|---|
| 154 | ! DECLARATION |
|---|
| 155 | ! ----------- |
|---|
| 156 | implicit none |
|---|
| 157 | |
|---|
| 158 | ! ARGUMENTS |
|---|
| 159 | ! --------- |
|---|
| 160 | real(dp), dimension(:,:,:), intent(in) :: q_PCM_in |
|---|
| 161 | integer(di), intent(in) :: iq |
|---|
| 162 | |
|---|
| 163 | ! LOCAL VARIABLES |
|---|
| 164 | ! --------------- |
|---|
| 165 | integer(di) :: l |
|---|
| 166 | |
|---|
| 167 | ! CODE |
|---|
| 168 | ! ---- |
|---|
| 169 | do l = 1,size(q_PCM_in,3) |
|---|
| 170 | call dyngrd2vect(nlon + 1,nlat,ngrid,q_PCM_in(:,:,l),q_PCM(:,l,iq)) |
|---|
| 171 | end do |
|---|
| 172 | |
|---|
| 173 | END SUBROUTINE set_q_PCM |
|---|
| 174 | !======================================================================= |
|---|
| 175 | |
|---|
| 176 | !======================================================================= |
|---|
| 177 | SUBROUTINE adapt_tracers2pressure(ps_avg_glob_old,ps_avg_glob,ps_ts,q_h2o_ts,q_co2_ts) |
|---|
| 178 | !----------------------------------------------------------------------- |
|---|
| 179 | ! NAME |
|---|
| 180 | ! adapt_tracers2pressure |
|---|
| 181 | ! |
|---|
| 182 | ! DESCRIPTION |
|---|
| 183 | ! Adapt the timeseries of tracers VMR to the new pressure. |
|---|
| 184 | ! |
|---|
| 185 | ! AUTHORS & DATE |
|---|
| 186 | ! JB Clement, 01/2026 |
|---|
| 187 | ! |
|---|
| 188 | ! NOTES |
|---|
| 189 | ! |
|---|
| 190 | !----------------------------------------------------------------------- |
|---|
| 191 | |
|---|
| 192 | ! DEPENDENCIES |
|---|
| 193 | ! ------------ |
|---|
| 194 | use geometry, only: ngrid, nday |
|---|
| 195 | use atmosphere, only: ap, bp |
|---|
| 196 | use display, only: print_msg |
|---|
| 197 | |
|---|
| 198 | ! DECLARATION |
|---|
| 199 | ! ----------- |
|---|
| 200 | implicit none |
|---|
| 201 | |
|---|
| 202 | ! ARGUMENTS |
|---|
| 203 | ! --------- |
|---|
| 204 | real(dp), intent(in) :: ps_avg_glob_old, ps_avg_glob |
|---|
| 205 | real(dp), dimension(ngrid,nday), intent(inout) :: ps_ts, q_h2o_ts, q_co2_ts |
|---|
| 206 | |
|---|
| 207 | ! LOCAL VARIABLES |
|---|
| 208 | ! --------------- |
|---|
| 209 | real(dp), dimension(ngrid,nday) :: plev1, plev2, plev1_old, plev2_old |
|---|
| 210 | |
|---|
| 211 | ! CODE |
|---|
| 212 | ! ---- |
|---|
| 213 | call print_msg("> Adapting the timeseries of tracers VMR to the new pressure") |
|---|
| 214 | |
|---|
| 215 | ! Build the first levels of pressure |
|---|
| 216 | plev1_old(:,:) = ap(1) + bp(1)*ps_ts(:,:) |
|---|
| 217 | plev2_old(:,:) = ap(2) + bp(2)*ps_ts(:,:) |
|---|
| 218 | |
|---|
| 219 | ! Evolve timeseries of surface pressure |
|---|
| 220 | ps_ts(:,:) = ps_ts(:,:)*ps_avg_glob/ps_avg_glob_old |
|---|
| 221 | |
|---|
| 222 | ! Build the first levels of pressure |
|---|
| 223 | plev1(:,:) = ap(1) + bp(1)*ps_ts(:,:) |
|---|
| 224 | plev2(:,:) = ap(2) + bp(2)*ps_ts(:,:) |
|---|
| 225 | |
|---|
| 226 | ! Adapt the timeseries of VMR tracers |
|---|
| 227 | ! H2O |
|---|
| 228 | q_h2o_ts(:,:) = q_h2o_ts(:,:)*(plev1_old(:,:) - plev2_old(:,:))/(plev1(:,:) - plev2(:,:)) |
|---|
| 229 | where (q_h2o_ts < 0._dp) q_h2o_ts = 0._dp |
|---|
| 230 | where (q_h2o_ts > 1._dp) q_h2o_ts = 1._dp |
|---|
| 231 | ! CO2 |
|---|
| 232 | q_co2_ts(:,:) = q_co2_ts(:,:)*(plev1_old(:,:) - plev2_old(:,:))/(plev1(:,:) - plev2(:,:)) |
|---|
| 233 | where (q_co2_ts < 0._dp) q_co2_ts = 0._dp |
|---|
| 234 | where (q_co2_ts > 1._dp) q_co2_ts = 1._dp |
|---|
| 235 | |
|---|
| 236 | END SUBROUTINE adapt_tracers2pressure |
|---|
| 237 | !======================================================================= |
|---|
| 238 | |
|---|
| 239 | !======================================================================= |
|---|
| 240 | SUBROUTINE build4PCM_tracers(ps4PCM,q4PCM) |
|---|
| 241 | !----------------------------------------------------------------------- |
|---|
| 242 | ! NAME |
|---|
| 243 | ! build4PCM_tracers |
|---|
| 244 | ! |
|---|
| 245 | ! DESCRIPTION |
|---|
| 246 | ! Reconstructs the tracer VMRs for the PCM. |
|---|
| 247 | ! |
|---|
| 248 | ! AUTHORS & DATE |
|---|
| 249 | ! JB Clement, 12/2025 |
|---|
| 250 | ! C. Metz, 01/2026 |
|---|
| 251 | ! |
|---|
| 252 | ! NOTES |
|---|
| 253 | ! |
|---|
| 254 | !----------------------------------------------------------------------- |
|---|
| 255 | |
|---|
| 256 | ! DEPENDENCIES |
|---|
| 257 | ! ------------ |
|---|
| 258 | use atmosphere, only: ap, bp, ps_PCM |
|---|
| 259 | use geometry, only: ngrid, nlayer |
|---|
| 260 | use display, only: print_msg |
|---|
| 261 | |
|---|
| 262 | ! DECLARATION |
|---|
| 263 | ! ----------- |
|---|
| 264 | implicit none |
|---|
| 265 | |
|---|
| 266 | ! ARGUMENTS |
|---|
| 267 | ! --------- |
|---|
| 268 | real(dp), dimension(:), intent(in) :: ps4PCM |
|---|
| 269 | real(dp), dimension(:,:,:), intent(out) :: q4PCM |
|---|
| 270 | |
|---|
| 271 | ! LOCAL VARIABLES |
|---|
| 272 | ! --------------- |
|---|
| 273 | integer(di) :: i, l, iq |
|---|
| 274 | real(dp) :: extra_mass |
|---|
| 275 | real(dp), dimension(:,:), allocatable :: dz_plev_PCM, dz_plev_new |
|---|
| 276 | |
|---|
| 277 | ! CODE |
|---|
| 278 | ! ---- |
|---|
| 279 | call print_msg('> Building tracers for the PCM') |
|---|
| 280 | allocate(dz_plev_PCM(ngrid,nlayer),dz_plev_new(ngrid,nlayer)) |
|---|
| 281 | |
|---|
| 282 | ! Compute thickness of pressure levels (= atmospheric layers) |
|---|
| 283 | do l = 1,nlayer |
|---|
| 284 | dz_plev_PCM(:,l) = ap(l) + bp(l)*ps_PCM(:) - ap(l + 1) - bp(l + 1)*ps_PCM(:) |
|---|
| 285 | dz_plev_new(:,l) = ap(l) + bp(l)*ps_PCM(:) - ap(l + 1) - bp(l + 1)*ps4PCM(:) |
|---|
| 286 | end do |
|---|
| 287 | |
|---|
| 288 | ! Reconstruct tracers |
|---|
| 289 | do iq = 1,nq |
|---|
| 290 | do l = 1,nlayer |
|---|
| 291 | do i = 1,ngrid |
|---|
| 292 | ! Standard scaling |
|---|
| 293 | q4PCM(i,l,iq) = q_PCM(i,l,iq)*dz_plev_PCM(i,l)/dz_plev_new(i,l) |
|---|
| 294 | |
|---|
| 295 | ! CO2 logic: account for the dry air mass change |
|---|
| 296 | if (trim(qnames(iq)) == "co2") q4PCM(i,l,iq) = q4PCM(i,l,iq) + (dz_plev_new(i,l) - dz_plev_PCM(i,l)) / dz_plev_new(i,l) |
|---|
| 297 | end do |
|---|
| 298 | end do |
|---|
| 299 | end do |
|---|
| 300 | |
|---|
| 301 | ! Mass conservation and clipping |
|---|
| 302 | do iq = 1,nq |
|---|
| 303 | do i = 1,ngrid |
|---|
| 304 | do l = 1,nlayer |
|---|
| 305 | ! Clipping negative values |
|---|
| 306 | if (q4PCM(i,l,iq) < 0._dp) q4PCM(i,l,iq) = 0._dp |
|---|
| 307 | |
|---|
| 308 | ! VMR limit check (ignore number density tracers) |
|---|
| 309 | if (q4PCM(i,l,iq) > 1._dp .and. & |
|---|
| 310 | (qnames(iq) /= "dust_number") .and. & |
|---|
| 311 | (qnames(iq) /= "ccn_number") .and. & |
|---|
| 312 | (qnames(iq) /= "stormdust_number") .and. & |
|---|
| 313 | (qnames(iq) /= "topdust_number")) then |
|---|
| 314 | |
|---|
| 315 | if (l < nlayer) then |
|---|
| 316 | extra_mass = (q4PCM(i,l,iq) - 1.)*dz_plev_new(i,l) ! extra_mass units: [VMR*pressure] |
|---|
| 317 | q4PCM(i,l,iq) = 1._dp |
|---|
| 318 | |
|---|
| 319 | ! Transfer extra mass to the layer above |
|---|
| 320 | ! Note: We divide by dz_plev_new(i,l + 1) to convert back to VMR |
|---|
| 321 | q4PCM(i,l + 1,iq) = q4PCM(i,l + 1,iq) + extra_mass/dz_plev_new(i,l + 1) |
|---|
| 322 | else |
|---|
| 323 | q4PCM(i,l,iq) = 1._dp |
|---|
| 324 | end if |
|---|
| 325 | end if |
|---|
| 326 | end do |
|---|
| 327 | end do |
|---|
| 328 | end do |
|---|
| 329 | |
|---|
| 330 | deallocate(dz_plev_PCM,dz_plev_new) |
|---|
| 331 | |
|---|
| 332 | END SUBROUTINE build4PCM_tracers |
|---|
| 333 | !======================================================================= |
|---|
| 334 | |
|---|
| 335 | END MODULE tracers |
|---|