| 1 | MODULE atmosphere |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! atmosphere |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! Contains global parameters used for the atmosphere. |
|---|
| 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(9), parameter :: z2sigdef_name = 'z2sig.def' |
|---|
| 27 | logical(k4), protected :: hybrid_alt_coord ! Flag fo hybrid coordinates |
|---|
| 28 | real(dp), dimension(:), allocatable, protected :: ap ! Hybrid (pressure contribution) coordinate at layer interfaces [Pa] |
|---|
| 29 | real(dp), dimension(:), allocatable, protected :: bp ! Hybrid (sigma contribution) coordinate at layer interfaces |
|---|
| 30 | real(dp), dimension(:), allocatable, protected :: ps_PCM ! Surface pressure in the "start.nc" at the beginning [Pa] |
|---|
| 31 | real(dp), dimension(:,:), allocatable, protected :: teta_PCM ! Potential temperature in the "start.nc" at the beginning [K] |
|---|
| 32 | real(dp), dimension(:,:), allocatable, protected :: u_PCM ! Zonal wind [m/s] |
|---|
| 33 | real(dp), dimension(:,:), allocatable, protected :: v_PCM ! Meridional wind [m/s] |
|---|
| 34 | real(dp) :: CO2cond_ps_PCM ! Coefficient to control the surface pressure change |
|---|
| 35 | |
|---|
| 36 | contains |
|---|
| 37 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 38 | |
|---|
| 39 | !======================================================================= |
|---|
| 40 | SUBROUTINE ini_atmosphere() |
|---|
| 41 | !----------------------------------------------------------------------- |
|---|
| 42 | ! NAME |
|---|
| 43 | ! ini_atmosphere |
|---|
| 44 | ! |
|---|
| 45 | ! DESCRIPTION |
|---|
| 46 | ! Initialize the parameters of module 'atmosphere'. |
|---|
| 47 | ! |
|---|
| 48 | ! AUTHORS & DATE |
|---|
| 49 | ! JB Clement, 12/2025 |
|---|
| 50 | ! |
|---|
| 51 | ! NOTES |
|---|
| 52 | ! |
|---|
| 53 | !----------------------------------------------------------------------- |
|---|
| 54 | |
|---|
| 55 | ! DEPENDENCIES |
|---|
| 56 | ! ------------ |
|---|
| 57 | use geometry, only: nlayer, ngrid |
|---|
| 58 | |
|---|
| 59 | ! DECLARATION |
|---|
| 60 | ! ----------- |
|---|
| 61 | implicit none |
|---|
| 62 | |
|---|
| 63 | ! CODE |
|---|
| 64 | ! ---- |
|---|
| 65 | if (.not. allocated(ap)) allocate(ap(nlayer + 1)) |
|---|
| 66 | if (.not. allocated(bp)) allocate(bp(nlayer + 1)) |
|---|
| 67 | if (.not. allocated(ps_PCM)) allocate(ps_PCM(ngrid)) |
|---|
| 68 | if (.not. allocated(teta_PCM)) allocate(teta_PCM(ngrid,nlayer)) |
|---|
| 69 | if (.not. allocated(u_PCM)) allocate(u_PCM(ngrid,nlayer)) |
|---|
| 70 | if (.not. allocated(v_PCM)) allocate(v_PCM(ngrid,nlayer)) |
|---|
| 71 | |
|---|
| 72 | END SUBROUTINE ini_atmosphere |
|---|
| 73 | !======================================================================= |
|---|
| 74 | |
|---|
| 75 | !======================================================================= |
|---|
| 76 | SUBROUTINE end_atmosphere() |
|---|
| 77 | !----------------------------------------------------------------------- |
|---|
| 78 | ! NAME |
|---|
| 79 | ! end_atmosphere |
|---|
| 80 | ! |
|---|
| 81 | ! DESCRIPTION |
|---|
| 82 | ! Deallocate atmosphere arrays. |
|---|
| 83 | ! |
|---|
| 84 | ! AUTHORS & DATE |
|---|
| 85 | ! JB Clement, 12/2025 |
|---|
| 86 | ! |
|---|
| 87 | ! NOTES |
|---|
| 88 | ! |
|---|
| 89 | !----------------------------------------------------------------------- |
|---|
| 90 | |
|---|
| 91 | ! DECLARATION |
|---|
| 92 | ! ----------- |
|---|
| 93 | implicit none |
|---|
| 94 | |
|---|
| 95 | ! CODE |
|---|
| 96 | ! ---- |
|---|
| 97 | if (allocated(ap)) deallocate(ap) |
|---|
| 98 | if (allocated(bp)) deallocate(bp) |
|---|
| 99 | if (allocated(ps_PCM)) deallocate(ps_PCM) |
|---|
| 100 | if (allocated(teta_PCM)) deallocate(teta_PCM) |
|---|
| 101 | if (allocated(u_PCM)) deallocate(u_PCM) |
|---|
| 102 | if (allocated(v_PCM)) deallocate(v_PCM) |
|---|
| 103 | |
|---|
| 104 | END SUBROUTINE end_atmosphere |
|---|
| 105 | !======================================================================= |
|---|
| 106 | |
|---|
| 107 | !======================================================================= |
|---|
| 108 | SUBROUTINE set_atmosphere_config(hybrid_alt_coord_in) |
|---|
| 109 | !----------------------------------------------------------------------- |
|---|
| 110 | ! NAME |
|---|
| 111 | ! set_atmosphere_config |
|---|
| 112 | ! |
|---|
| 113 | ! DESCRIPTION |
|---|
| 114 | ! Setter for 'atmosphere' configuration parameters. |
|---|
| 115 | ! |
|---|
| 116 | ! AUTHORS & DATE |
|---|
| 117 | ! JB Clement, 02/2026 |
|---|
| 118 | ! |
|---|
| 119 | ! NOTES |
|---|
| 120 | ! |
|---|
| 121 | !----------------------------------------------------------------------- |
|---|
| 122 | |
|---|
| 123 | ! DEPENDENCIES |
|---|
| 124 | ! ------------ |
|---|
| 125 | use utility, only: bool2str |
|---|
| 126 | use display, only: print_msg |
|---|
| 127 | |
|---|
| 128 | ! DECLARATION |
|---|
| 129 | ! ----------- |
|---|
| 130 | implicit none |
|---|
| 131 | |
|---|
| 132 | ! ARGUMENTS |
|---|
| 133 | ! --------- |
|---|
| 134 | logical(k4), intent(in) :: hybrid_alt_coord_in |
|---|
| 135 | |
|---|
| 136 | ! CODE |
|---|
| 137 | ! ---- |
|---|
| 138 | hybrid_alt_coord = hybrid_alt_coord_in |
|---|
| 139 | call print_msg('hybrid_alt_coord = '//bool2str(hybrid_alt_coord)) |
|---|
| 140 | |
|---|
| 141 | END SUBROUTINE set_atmosphere_config |
|---|
| 142 | !======================================================================= |
|---|
| 143 | |
|---|
| 144 | !======================================================================= |
|---|
| 145 | SUBROUTINE set_ap(ap_in) |
|---|
| 146 | !----------------------------------------------------------------------- |
|---|
| 147 | ! NAME |
|---|
| 148 | ! set_ap |
|---|
| 149 | ! |
|---|
| 150 | ! DESCRIPTION |
|---|
| 151 | ! Setter for 'ap'. |
|---|
| 152 | ! |
|---|
| 153 | ! AUTHORS & DATE |
|---|
| 154 | ! JB Clement, 12/2025 |
|---|
| 155 | ! |
|---|
| 156 | ! NOTES |
|---|
| 157 | ! |
|---|
| 158 | !----------------------------------------------------------------------- |
|---|
| 159 | |
|---|
| 160 | ! DECLARATION |
|---|
| 161 | ! ----------- |
|---|
| 162 | implicit none |
|---|
| 163 | |
|---|
| 164 | ! ARGUMENTS |
|---|
| 165 | ! --------- |
|---|
| 166 | real(dp), dimension(:), intent(in) :: ap_in |
|---|
| 167 | |
|---|
| 168 | ! CODE |
|---|
| 169 | ! ---- |
|---|
| 170 | ap(:) = ap_in(:) |
|---|
| 171 | |
|---|
| 172 | END SUBROUTINE set_ap |
|---|
| 173 | !======================================================================= |
|---|
| 174 | |
|---|
| 175 | !======================================================================= |
|---|
| 176 | SUBROUTINE set_bp(bp_in) |
|---|
| 177 | !----------------------------------------------------------------------- |
|---|
| 178 | ! NAME |
|---|
| 179 | ! set_bp |
|---|
| 180 | ! |
|---|
| 181 | ! DESCRIPTION |
|---|
| 182 | ! Setter for 'bp'. |
|---|
| 183 | ! |
|---|
| 184 | ! AUTHORS & DATE |
|---|
| 185 | ! JB Clement, 12/2025 |
|---|
| 186 | ! |
|---|
| 187 | ! NOTES |
|---|
| 188 | ! |
|---|
| 189 | !----------------------------------------------------------------------- |
|---|
| 190 | |
|---|
| 191 | ! DECLARATION |
|---|
| 192 | ! ----------- |
|---|
| 193 | implicit none |
|---|
| 194 | |
|---|
| 195 | ! ARGUMENTS |
|---|
| 196 | ! --------- |
|---|
| 197 | real(dp), dimension(:), intent(in) :: bp_in |
|---|
| 198 | |
|---|
| 199 | ! CODE |
|---|
| 200 | ! ---- |
|---|
| 201 | bp(:) = bp_in(:) |
|---|
| 202 | |
|---|
| 203 | END SUBROUTINE set_bp |
|---|
| 204 | !======================================================================= |
|---|
| 205 | |
|---|
| 206 | !======================================================================= |
|---|
| 207 | SUBROUTINE set_ps_PCM(ps_PCM_in) |
|---|
| 208 | !----------------------------------------------------------------------- |
|---|
| 209 | ! NAME |
|---|
| 210 | ! set_ps_PCM |
|---|
| 211 | ! |
|---|
| 212 | ! DESCRIPTION |
|---|
| 213 | ! Setter for 'ps_PCM'. |
|---|
| 214 | ! |
|---|
| 215 | ! AUTHORS & DATE |
|---|
| 216 | ! JB Clement, 12/2025 |
|---|
| 217 | ! |
|---|
| 218 | ! NOTES |
|---|
| 219 | ! |
|---|
| 220 | !----------------------------------------------------------------------- |
|---|
| 221 | |
|---|
| 222 | ! DEPENDENCIES |
|---|
| 223 | ! ------------ |
|---|
| 224 | use geometry, only: dyngrd2vect, nlon, nlat, ngrid |
|---|
| 225 | |
|---|
| 226 | ! DECLARATION |
|---|
| 227 | ! ----------- |
|---|
| 228 | implicit none |
|---|
| 229 | |
|---|
| 230 | ! ARGUMENTS |
|---|
| 231 | ! --------- |
|---|
| 232 | real(dp), dimension(:,:), intent(in) :: ps_PCM_in |
|---|
| 233 | |
|---|
| 234 | ! CODE |
|---|
| 235 | ! ---- |
|---|
| 236 | call dyngrd2vect(nlon + 1,nlat,ngrid,ps_PCM_in,ps_PCM) |
|---|
| 237 | |
|---|
| 238 | END SUBROUTINE set_ps_PCM |
|---|
| 239 | !======================================================================= |
|---|
| 240 | |
|---|
| 241 | !======================================================================= |
|---|
| 242 | SUBROUTINE set_teta_PCM(teta_PCM_in) |
|---|
| 243 | !----------------------------------------------------------------------- |
|---|
| 244 | ! NAME |
|---|
| 245 | ! set_teta_PCM |
|---|
| 246 | ! |
|---|
| 247 | ! DESCRIPTION |
|---|
| 248 | ! Setter for 'teta_PCM'. |
|---|
| 249 | ! |
|---|
| 250 | ! AUTHORS & DATE |
|---|
| 251 | ! JB Clement, 01/2026 |
|---|
| 252 | ! |
|---|
| 253 | ! NOTES |
|---|
| 254 | ! |
|---|
| 255 | !----------------------------------------------------------------------- |
|---|
| 256 | |
|---|
| 257 | ! DEPENDENCIES |
|---|
| 258 | ! ------------ |
|---|
| 259 | use geometry, only: dyngrd2vect, nlon, nlat, ngrid, nlayer |
|---|
| 260 | |
|---|
| 261 | ! DECLARATION |
|---|
| 262 | ! ----------- |
|---|
| 263 | implicit none |
|---|
| 264 | |
|---|
| 265 | ! ARGUMENTS |
|---|
| 266 | ! --------- |
|---|
| 267 | real(dp), dimension(:,:,:), intent(in) :: teta_PCM_in |
|---|
| 268 | |
|---|
| 269 | ! LOCAL VARIABLES |
|---|
| 270 | ! --------------- |
|---|
| 271 | integer(di) :: l |
|---|
| 272 | |
|---|
| 273 | ! CODE |
|---|
| 274 | ! ---- |
|---|
| 275 | do l = 1,nlayer |
|---|
| 276 | call dyngrd2vect(nlon + 1,nlat,ngrid,teta_PCM_in(:,:,l),teta_PCM(:,l)) |
|---|
| 277 | end do |
|---|
| 278 | |
|---|
| 279 | END SUBROUTINE set_teta_PCM |
|---|
| 280 | !======================================================================= |
|---|
| 281 | |
|---|
| 282 | !======================================================================= |
|---|
| 283 | SUBROUTINE set_u_PCM(u_PCM_in) |
|---|
| 284 | !----------------------------------------------------------------------- |
|---|
| 285 | ! NAME |
|---|
| 286 | ! set_u_PCM |
|---|
| 287 | ! |
|---|
| 288 | ! DESCRIPTION |
|---|
| 289 | ! Setter for 'u_PCM'. |
|---|
| 290 | ! |
|---|
| 291 | ! AUTHORS & DATE |
|---|
| 292 | ! JB Clement, 01/2026 |
|---|
| 293 | ! |
|---|
| 294 | ! NOTES |
|---|
| 295 | ! |
|---|
| 296 | !----------------------------------------------------------------------- |
|---|
| 297 | |
|---|
| 298 | ! DEPENDENCIES |
|---|
| 299 | ! ------------ |
|---|
| 300 | use geometry, only: dyngrd2vect, nlon, nlat, ngrid, nlayer |
|---|
| 301 | |
|---|
| 302 | ! DECLARATION |
|---|
| 303 | ! ----------- |
|---|
| 304 | implicit none |
|---|
| 305 | |
|---|
| 306 | ! ARGUMENTS |
|---|
| 307 | ! --------- |
|---|
| 308 | real(dp), dimension(:,:,:), intent(in) :: u_PCM_in |
|---|
| 309 | |
|---|
| 310 | ! LOCAL VARIABLES |
|---|
| 311 | ! --------------- |
|---|
| 312 | integer(di) :: l |
|---|
| 313 | |
|---|
| 314 | ! CODE |
|---|
| 315 | ! ---- |
|---|
| 316 | do l = 1,nlayer |
|---|
| 317 | call dyngrd2vect(nlon + 1,nlat,ngrid,u_PCM_in(:,:,l),u_PCM(:,l)) |
|---|
| 318 | end do |
|---|
| 319 | |
|---|
| 320 | END SUBROUTINE set_u_PCM |
|---|
| 321 | !======================================================================= |
|---|
| 322 | |
|---|
| 323 | !======================================================================= |
|---|
| 324 | SUBROUTINE set_v_PCM(v_PCM_in) |
|---|
| 325 | !----------------------------------------------------------------------- |
|---|
| 326 | ! NAME |
|---|
| 327 | ! set_v_PCM |
|---|
| 328 | ! |
|---|
| 329 | ! DESCRIPTION |
|---|
| 330 | ! Setter for 'v_PCM'. |
|---|
| 331 | ! |
|---|
| 332 | ! AUTHORS & DATE |
|---|
| 333 | ! JB Clement, 01/2026 |
|---|
| 334 | ! |
|---|
| 335 | ! NOTES |
|---|
| 336 | ! |
|---|
| 337 | !----------------------------------------------------------------------- |
|---|
| 338 | |
|---|
| 339 | ! DEPENDENCIES |
|---|
| 340 | ! ------------ |
|---|
| 341 | use geometry, only: dyngrd2vect, nlon, nlat, ngrid, nlayer |
|---|
| 342 | |
|---|
| 343 | ! DECLARATION |
|---|
| 344 | ! ----------- |
|---|
| 345 | implicit none |
|---|
| 346 | |
|---|
| 347 | ! ARGUMENTS |
|---|
| 348 | ! --------- |
|---|
| 349 | real(dp), dimension(:,:,:), intent(in) :: v_PCM_in |
|---|
| 350 | |
|---|
| 351 | ! LOCAL VARIABLES |
|---|
| 352 | ! --------------- |
|---|
| 353 | integer(di) :: l |
|---|
| 354 | |
|---|
| 355 | ! CODE |
|---|
| 356 | ! ---- |
|---|
| 357 | do l = 1,nlayer |
|---|
| 358 | call dyngrd2vect(nlon + 1,nlat,ngrid,v_PCM_in(:,:,l),v_PCM(:,l)) |
|---|
| 359 | end do |
|---|
| 360 | |
|---|
| 361 | END SUBROUTINE set_v_PCM |
|---|
| 362 | !======================================================================= |
|---|
| 363 | |
|---|
| 364 | !======================================================================= |
|---|
| 365 | SUBROUTINE compute_alt_coord(pa,preff,ap,bp,aps,bps) |
|---|
| 366 | !----------------------------------------------------------------------- |
|---|
| 367 | ! NAME |
|---|
| 368 | ! compute_alt_coord |
|---|
| 369 | ! |
|---|
| 370 | ! DESCRIPTION |
|---|
| 371 | ! Compute values of the coefficients to define the altitude |
|---|
| 372 | ! coordinates. |
|---|
| 373 | ! |
|---|
| 374 | ! AUTHORS & DATE |
|---|
| 375 | ! JB Clement, 02/2026 |
|---|
| 376 | ! |
|---|
| 377 | ! NOTES |
|---|
| 378 | ! Taken from 'disvert_noterre' written by F. Forget, Y. Wanherdrick |
|---|
| 379 | ! and P. Levan in the dynamics. |
|---|
| 380 | ! |
|---|
| 381 | ! 'pa' and 'preff' are read from "start.nc"/"start1D.txt". |
|---|
| 382 | ! |
|---|
| 383 | ! WARNING: to compute the values at mid_layer, there is an arbitrary |
|---|
| 384 | ! decision here to put the middle of the layer half-way (pressure-wise) |
|---|
| 385 | ! between boundaries. |
|---|
| 386 | ! In addition, for the last layer (the top of which is at zero pressure) |
|---|
| 387 | ! P(llm), we choose to put it with the same interval (in log(pressure)) |
|---|
| 388 | ! from P(llm-1), than there is between between P(llm-1) and P(llm-2) |
|---|
| 389 | ! This choice is quite specific and must be compliant with the PCM! |
|---|
| 390 | !----------------------------------------------------------------------- |
|---|
| 391 | |
|---|
| 392 | ! DEPENDENCIES |
|---|
| 393 | ! ------------ |
|---|
| 394 | use geometry, only: nlayer |
|---|
| 395 | use display, only: print_msg |
|---|
| 396 | use stoppage, only: stop_clean |
|---|
| 397 | use utility, only: real2str |
|---|
| 398 | |
|---|
| 399 | ! DECLARATION |
|---|
| 400 | ! ----------- |
|---|
| 401 | implicit none |
|---|
| 402 | |
|---|
| 403 | ! ARGUMENTS |
|---|
| 404 | ! --------- |
|---|
| 405 | real(dp) :: pa, preff |
|---|
| 406 | real(dp), dimension(nlayer + 1), intent(out) :: ap, bp |
|---|
| 407 | real(dp), dimension(nlayer), intent(out), optional :: aps, bps |
|---|
| 408 | |
|---|
| 409 | ! LOCAL VAIRABLES |
|---|
| 410 | ! --------------- |
|---|
| 411 | integer(di) :: i |
|---|
| 412 | real(dp) :: scaleheight, newsig |
|---|
| 413 | real(dp), dimension(nlayer + 1) :: sig |
|---|
| 414 | |
|---|
| 415 | ! CODE |
|---|
| 416 | ! ---- |
|---|
| 417 | ! Read vertical grid definition |
|---|
| 418 | call read_z2sig(scaleheight,sig) |
|---|
| 419 | call print_msg('Pa [Pa] = '//real2str(pa)) |
|---|
| 420 | call print_msg('Preff [Pa] = '//real2str(preff)) |
|---|
| 421 | |
|---|
| 422 | ! Compute ap and bp coefficients |
|---|
| 423 | if (hybrid_alt_coord) then |
|---|
| 424 | call print_msg('> Defining hybrid altitude coordinates') |
|---|
| 425 | do i = 1,nlayer |
|---|
| 426 | call compute_hybrid_sig(sig(i),pa,preff,newsig) |
|---|
| 427 | bp(i) = exp(1._dp - 1._dp/(newsig**2)) |
|---|
| 428 | ap(i) = pa*(newsig - bp(i)) |
|---|
| 429 | end do |
|---|
| 430 | ap(nlayer + 1) = 0._dp |
|---|
| 431 | bp(nlayer + 1) = 0._dp |
|---|
| 432 | else |
|---|
| 433 | call print_msg('> Defining sigma altitude coordinates') |
|---|
| 434 | ap(:) = 0._dp |
|---|
| 435 | bp(1:nlayer) = sig(1:nlayer) |
|---|
| 436 | bp(nlayer + 1) = 0._dp |
|---|
| 437 | end if |
|---|
| 438 | |
|---|
| 439 | ! Compute aps and bps coefficients (mid-layer) if asked |
|---|
| 440 | if (present(aps) .neqv. present(bps)) call stop_clean(__FILE__,__LINE__,"'aps' and 'bps' must be computed together!",1) |
|---|
| 441 | if (.not. present(aps)) return |
|---|
| 442 | |
|---|
| 443 | do i = 1,nlayer - 1 |
|---|
| 444 | aps(i) = 0.5_dp*(ap(i) + ap(i + 1)) |
|---|
| 445 | bps(i) = 0.5_dp*(bp(i) + bp(i + 1)) |
|---|
| 446 | end do |
|---|
| 447 | if (hybrid_alt_coord) then |
|---|
| 448 | aps(nlayer) = aps(nlayer - 1)**2/aps(nlayer - 2) |
|---|
| 449 | bps(nlayer) = 0.5_dp*(bp(nlayer) + bp(nlayer + 1)) |
|---|
| 450 | else |
|---|
| 451 | aps(nlayer) = 0._dp |
|---|
| 452 | bps(nlayer) = bps(nlayer - 1)**2/bps(nlayer - 2) |
|---|
| 453 | end if |
|---|
| 454 | |
|---|
| 455 | END SUBROUTINE compute_alt_coord |
|---|
| 456 | !======================================================================= |
|---|
| 457 | |
|---|
| 458 | !======================================================================= |
|---|
| 459 | SUBROUTINE compute_hybrid_sig(sig,pa,preff,newsig) |
|---|
| 460 | !----------------------------------------------------------------------- |
|---|
| 461 | ! NAME |
|---|
| 462 | ! compute_hybrid_sig |
|---|
| 463 | ! |
|---|
| 464 | ! DESCRIPTION |
|---|
| 465 | ! Compute modified sigma value to keep vertical coordinates described |
|---|
| 466 | ! in "z2sig.def" when defining hybrid coordinates. |
|---|
| 467 | ! |
|---|
| 468 | ! AUTHORS & DATE |
|---|
| 469 | ! JB Clement, 01/2026 |
|---|
| 470 | ! |
|---|
| 471 | ! NOTES |
|---|
| 472 | ! Taken from 'sig_hybrid' written by F. Forget (2002) in the |
|---|
| 473 | ! dynamics. |
|---|
| 474 | ! |
|---|
| 475 | ! Knowing sig (the "sigma" levels where we want to put the layers), |
|---|
| 476 | ! it solves 'newsig' such that: |
|---|
| 477 | ! (1 - pa/preff)*exp(1 - 1/newsig^2) + (pa/preff)*newsig = sig |
|---|
| 478 | ! which cannot be done analyticaly. |
|---|
| 479 | ! => needs for an iterative procedure |
|---|
| 480 | ! |
|---|
| 481 | ! Information: when exp(1 - 1/x**2) becomes << x |
|---|
| 482 | ! x exp(1 - 1/x**2)/x |
|---|
| 483 | ! 1 1 |
|---|
| 484 | ! 0.68 0.5 |
|---|
| 485 | ! 0.5 1.E-1 |
|---|
| 486 | ! 0.391 1.E-2 |
|---|
| 487 | ! 0.333 1.E-3 |
|---|
| 488 | ! 0.295 1.E-4 |
|---|
| 489 | ! 0.269 1.E-5 |
|---|
| 490 | ! 0.248 1.E-6 |
|---|
| 491 | ! => one can thus use newsig = sig*preff/pa if sig*preff/pa < 0.25 |
|---|
| 492 | !----------------------------------------------------------------------- |
|---|
| 493 | |
|---|
| 494 | ! DEPENDENCIES |
|---|
| 495 | ! ------------ |
|---|
| 496 | use stoppage, only: stop_clean |
|---|
| 497 | |
|---|
| 498 | ! DECLARATION |
|---|
| 499 | ! ----------- |
|---|
| 500 | implicit none |
|---|
| 501 | |
|---|
| 502 | ! ARGUMENTS |
|---|
| 503 | ! --------- |
|---|
| 504 | real(dp), intent(in) :: sig, pa, preff |
|---|
| 505 | real(dp), intent(out) :: newsig |
|---|
| 506 | |
|---|
| 507 | ! LOCAL VARIABLES |
|---|
| 508 | ! --------------- |
|---|
| 509 | real(dp) :: x1, x2, F |
|---|
| 510 | integer(di) :: j |
|---|
| 511 | |
|---|
| 512 | ! CODE |
|---|
| 513 | ! ---- |
|---|
| 514 | newsig = sig |
|---|
| 515 | x1 = 0._dp |
|---|
| 516 | x2 = 1._dp |
|---|
| 517 | |
|---|
| 518 | if (sig >= 1.) then |
|---|
| 519 | newsig = sig |
|---|
| 520 | else if (sig*preff/pa >= 0.25_dp) then |
|---|
| 521 | do j = 1,9999 ! Overkill maximum number of iterations |
|---|
| 522 | F = ((1._dp - pa/preff)*exp(1._dp - 1._dp/newsig**2) + (pa/preff)*newsig)/sig |
|---|
| 523 | |
|---|
| 524 | if (F > 1.) then |
|---|
| 525 | x2 = newsig |
|---|
| 526 | newsig = 0.5_dp*(x1 + newsig) |
|---|
| 527 | else |
|---|
| 528 | x1 = newsig |
|---|
| 529 | newsig = 0.5_dp*(x2 + newsig) |
|---|
| 530 | end if |
|---|
| 531 | ! Convergence is assumed if sig is within 0.01m (in pseudo-altitude) of the target value |
|---|
| 532 | if (abs(10._dp * log(F)) < 1.e-5_dp) exit |
|---|
| 533 | end do |
|---|
| 534 | else !if (sig*preff/pa <= 0.25_dp) then |
|---|
| 535 | newsig = sig*preff/pa |
|---|
| 536 | end if |
|---|
| 537 | |
|---|
| 538 | END SUBROUTINE compute_hybrid_sig |
|---|
| 539 | !======================================================================= |
|---|
| 540 | |
|---|
| 541 | !======================================================================= |
|---|
| 542 | SUBROUTINE read_z2sig(scaleheight,sig) |
|---|
| 543 | !----------------------------------------------------------------------- |
|---|
| 544 | ! NAME |
|---|
| 545 | ! read_z2sig |
|---|
| 546 | ! |
|---|
| 547 | ! DESCRIPTION |
|---|
| 548 | ! Read vertical grid definition from z2sig.def and compute sigma |
|---|
| 549 | ! levels. |
|---|
| 550 | ! |
|---|
| 551 | ! AUTHORS & DATE |
|---|
| 552 | ! JB Clement, 01/2026 |
|---|
| 553 | ! |
|---|
| 554 | ! NOTES |
|---|
| 555 | ! Adapted from 'disvert_noterre' written by F. Forget, Y. |
|---|
| 556 | ! Wanherdrick and P. Levan in the dynamics. |
|---|
| 557 | !----------------------------------------------------------------------- |
|---|
| 558 | |
|---|
| 559 | ! DEPENDENCIES |
|---|
| 560 | ! ------------ |
|---|
| 561 | use geometry, only: nlayer |
|---|
| 562 | use stoppage, only: stop_clean |
|---|
| 563 | use display, only: print_msg |
|---|
| 564 | use utility, only: real2str |
|---|
| 565 | |
|---|
| 566 | ! DECLARATION |
|---|
| 567 | ! ----------- |
|---|
| 568 | implicit none |
|---|
| 569 | |
|---|
| 570 | ! ARGUMENTS |
|---|
| 571 | ! --------- |
|---|
| 572 | real(dp), intent(out) :: scaleheight |
|---|
| 573 | real(dp), dimension(nlayer + 1), intent(out) :: sig |
|---|
| 574 | |
|---|
| 575 | ! LOCAL VARIABLES |
|---|
| 576 | ! --------------- |
|---|
| 577 | real(dp), dimension(nlayer) :: zsig |
|---|
| 578 | integer(di) :: l, ierr, funit |
|---|
| 579 | logical(k4) :: here |
|---|
| 580 | |
|---|
| 581 | ! CODE |
|---|
| 582 | ! ---- |
|---|
| 583 | ! Open "z2sig.def" |
|---|
| 584 | inquire(file = z2sigdef_name,exist = here) |
|---|
| 585 | if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//z2sigdef_name//'"!',1) |
|---|
| 586 | call print_msg('> Reading "'//z2sigdef_name//'"') |
|---|
| 587 | open(newunit = funit,file = z2sigdef_name,status = "old",action = "read",iostat = ierr) |
|---|
| 588 | if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening "'//z2sigdef_name//'"',ierr) |
|---|
| 589 | |
|---|
| 590 | ! Read scale height |
|---|
| 591 | read(funit,*,iostat = ierr) scaleheight |
|---|
| 592 | if (ierr /= 0) call stop_clean(__FILE__,__LINE__,"error reading 'scaleheight'!",1) |
|---|
| 593 | |
|---|
| 594 | ! Read pseudo-altitudes |
|---|
| 595 | do l = 1,nlayer |
|---|
| 596 | read(funit,*,iostat = ierr) zsig(l) |
|---|
| 597 | if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'"'//z2sigdef_name//'" too short for the number of altitude levels!',1) |
|---|
| 598 | end do |
|---|
| 599 | |
|---|
| 600 | ! Check for extra lines in the file (warning only) |
|---|
| 601 | read(funit,*,iostat = ierr) |
|---|
| 602 | if (ierr == 0) call print_msg('Warning: "'//z2sigdef_name//'" has more lines than expected!') |
|---|
| 603 | |
|---|
| 604 | ! Compute sigma values |
|---|
| 605 | sig(1) = 1. |
|---|
| 606 | do l = 2,nlayer |
|---|
| 607 | sig(l) = 0.5_dp*(exp(-zsig(l)/scaleheight) + exp(-zsig(l - 1)/scaleheight)) |
|---|
| 608 | end do |
|---|
| 609 | sig(nlayer + 1) = 0._dp |
|---|
| 610 | |
|---|
| 611 | ! Close |
|---|
| 612 | close(funit) |
|---|
| 613 | call print_msg('Scale height [km] = '//real2str(scaleheight)) |
|---|
| 614 | |
|---|
| 615 | END SUBROUTINE read_z2sig |
|---|
| 616 | !======================================================================= |
|---|
| 617 | |
|---|
| 618 | !======================================================================= |
|---|
| 619 | SUBROUTINE build4PCM_atmosphere(ps_avg,ps_dev,ps_avg_glob,ps_avg_glob_ini,ps4PCM,pa4PCM,preff4PCM,teta4PCM,air_mass4PCM) |
|---|
| 620 | !----------------------------------------------------------------------- |
|---|
| 621 | ! NAME |
|---|
| 622 | ! build4PCM_atmosphere |
|---|
| 623 | ! |
|---|
| 624 | ! DESCRIPTION |
|---|
| 625 | ! Reconstructs the pressure, potential temperature and tha air mass |
|---|
| 626 | ! for the PCM. |
|---|
| 627 | ! |
|---|
| 628 | ! AUTHORS & DATE |
|---|
| 629 | ! JB Clement, 01/2026 |
|---|
| 630 | ! |
|---|
| 631 | ! NOTES |
|---|
| 632 | ! Air mass computation taken from 'massdair' and 'pression' written |
|---|
| 633 | ! by P. Le Van and F. Hourdin in the dynamics. |
|---|
| 634 | !----------------------------------------------------------------------- |
|---|
| 635 | |
|---|
| 636 | ! DEPENDENCIES |
|---|
| 637 | ! ------------ |
|---|
| 638 | use geometry, only: ngrid, nlayer, cell_area, vect2dyngrd |
|---|
| 639 | use physics, only: g, rcp |
|---|
| 640 | use display, only: print_msg |
|---|
| 641 | |
|---|
| 642 | ! DECLARATION |
|---|
| 643 | ! ----------- |
|---|
| 644 | implicit none |
|---|
| 645 | |
|---|
| 646 | ! ARGUMENTS |
|---|
| 647 | ! --------- |
|---|
| 648 | real(dp), dimension(:), intent(in) :: ps_avg, ps_dev |
|---|
| 649 | real(dp), intent(in) :: ps_avg_glob, ps_avg_glob_ini |
|---|
| 650 | real(dp), dimension(:,:), intent(out) :: air_mass4PCM, teta4PCM |
|---|
| 651 | real(dp), dimension(:), intent(out) :: ps4PCM |
|---|
| 652 | real(dp), intent(out) :: pa4PCM, preff4PCM |
|---|
| 653 | |
|---|
| 654 | ! LOCAL VARIABLES |
|---|
| 655 | ! --------------- |
|---|
| 656 | real(dp), dimension(ngrid,nlayer + 1) :: p ! 3D pressure field |
|---|
| 657 | integer(di) :: l |
|---|
| 658 | |
|---|
| 659 | ! CODE |
|---|
| 660 | ! ---- |
|---|
| 661 | call print_msg('> Building pressure for the PCM') |
|---|
| 662 | ! The pressure deviation is rescaled to avoid disproportionate oscillations in case of huge change of average pressure during the PEM run |
|---|
| 663 | ps4PCM(:) = ps_avg(:) + ps_dev(:)*ps_avg_glob/ps_avg_glob_ini |
|---|
| 664 | pa4PCM = ps_avg_glob_ini/30. ! For now the altitude grid is not changed |
|---|
| 665 | preff4PCM = ps_avg_glob_ini ! For now the altitude grid is not changed |
|---|
| 666 | |
|---|
| 667 | call print_msg('> Building potential temperature for the PCM') |
|---|
| 668 | do l = 1,nlayer |
|---|
| 669 | ! Correction on teta due to surface pressure changes |
|---|
| 670 | teta4PCM(:,l) = teta_PCM(:,l)*ps4PCM(:)**rcp |
|---|
| 671 | end do |
|---|
| 672 | |
|---|
| 673 | call print_msg('> Building air mass for the PCM') |
|---|
| 674 | ! Compute atmospheric pressure |
|---|
| 675 | do l = 1,nlayer + 1 |
|---|
| 676 | p(:,l) = ap(l) + bp(l)*ps4PCM(:) |
|---|
| 677 | end do |
|---|
| 678 | |
|---|
| 679 | ! Compute air mass |
|---|
| 680 | do l = 1,nlayer |
|---|
| 681 | air_mass4PCM(:,l) = cell_area(:)*(p(:,l) - p(:,l + 1))/g |
|---|
| 682 | end do |
|---|
| 683 | |
|---|
| 684 | END SUBROUTINE build4PCM_atmosphere |
|---|
| 685 | !======================================================================= |
|---|
| 686 | |
|---|
| 687 | !======================================================================= |
|---|
| 688 | SUBROUTINE evolve_pressure(d_co2ice,delta_co2_ads,do_sorption,ps_avg_glob_old,ps_avg_glob,ps_avg) |
|---|
| 689 | !----------------------------------------------------------------------- |
|---|
| 690 | ! NAME |
|---|
| 691 | ! evolve_pressure |
|---|
| 692 | ! |
|---|
| 693 | ! DESCRIPTION |
|---|
| 694 | ! Evolve pressure according to tendency. |
|---|
| 695 | ! |
|---|
| 696 | ! AUTHORS & DATE |
|---|
| 697 | ! JB Clement, 01/2026 |
|---|
| 698 | ! |
|---|
| 699 | ! NOTES |
|---|
| 700 | ! |
|---|
| 701 | !----------------------------------------------------------------------- |
|---|
| 702 | |
|---|
| 703 | ! DEPENDENCIES |
|---|
| 704 | ! ------------ |
|---|
| 705 | use geometry, only: ngrid, nslope, cell_area, total_surface |
|---|
| 706 | use slopes, only: subslope_dist, def_slope_mean |
|---|
| 707 | use physics, only: g |
|---|
| 708 | use maths, only: pi |
|---|
| 709 | use display, only: print_msg |
|---|
| 710 | use utility, only: real2str |
|---|
| 711 | |
|---|
| 712 | ! DECLARATION |
|---|
| 713 | ! ----------- |
|---|
| 714 | implicit none |
|---|
| 715 | |
|---|
| 716 | ! ARGUMENTS |
|---|
| 717 | ! --------- |
|---|
| 718 | real(dp), dimension(:,:), intent(in) :: d_co2ice |
|---|
| 719 | real(dp), dimension(:), intent(in) :: delta_co2_ads |
|---|
| 720 | logical(k4), intent(in) :: do_sorption |
|---|
| 721 | real(dp), intent(inout) :: ps_avg_glob_old, ps_avg_glob |
|---|
| 722 | real(dp), dimension(:), intent(inout) :: ps_avg |
|---|
| 723 | |
|---|
| 724 | ! LOCAL VARIABLES |
|---|
| 725 | ! --------------- |
|---|
| 726 | integer(di) :: i, islope |
|---|
| 727 | |
|---|
| 728 | ! CODE |
|---|
| 729 | ! ---- |
|---|
| 730 | call print_msg("> Evolving the surface pressure") |
|---|
| 731 | ps_avg_glob_old = ps_avg_glob |
|---|
| 732 | do i = 1,ngrid |
|---|
| 733 | do islope = 1,nslope |
|---|
| 734 | ps_avg_glob = ps_avg_glob - CO2cond_ps_PCM*g*cell_area(i)*d_co2ice(i,islope)*subslope_dist(i,islope)/cos(pi*def_slope_mean(islope)/180.)/total_surface |
|---|
| 735 | end do |
|---|
| 736 | if (do_sorption) ps_avg_glob = ps_avg_glob - g*cell_area(i)*delta_co2_ads(i)/total_surface |
|---|
| 737 | end do |
|---|
| 738 | call print_msg('Global average pressure (old|new) [Pa] = '//real2str(ps_avg_glob_old)//' | '//real2str(ps_avg_glob)) |
|---|
| 739 | ps_avg(:) = ps_avg(:)*ps_avg_glob/ps_avg_glob_old |
|---|
| 740 | |
|---|
| 741 | END SUBROUTINE evolve_pressure |
|---|
| 742 | !======================================================================= |
|---|
| 743 | |
|---|
| 744 | END MODULE atmosphere |
|---|