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