| [4065] | 1 | MODULE surface |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! surface |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! Contains global parameters used for the surface. |
|---|
| 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 | real(dp), dimension(:), allocatable, protected :: albedodat_PCM ! Albedo of bare ground |
|---|
| 27 | real(dp), dimension(:,:), allocatable, protected :: albedo_PCM ! Surface albedo_PCM |
|---|
| 28 | real(dp), dimension(:,:), allocatable, protected :: emissivity_PCM ! Thermal IR surface emissivity_PCM |
|---|
| 29 | |
|---|
| 30 | contains |
|---|
| 31 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 32 | |
|---|
| 33 | !======================================================================= |
|---|
| 34 | SUBROUTINE ini_surface() |
|---|
| 35 | !----------------------------------------------------------------------- |
|---|
| 36 | ! NAME |
|---|
| 37 | ! ini_surface |
|---|
| 38 | ! |
|---|
| 39 | ! DESCRIPTION |
|---|
| 40 | ! Initialize the parameters of module 'surface'. |
|---|
| 41 | ! |
|---|
| 42 | ! AUTHORS & DATE |
|---|
| 43 | ! JB Clement, 12/2025 |
|---|
| 44 | ! |
|---|
| 45 | ! NOTES |
|---|
| 46 | ! |
|---|
| 47 | !----------------------------------------------------------------------- |
|---|
| 48 | |
|---|
| 49 | ! DEPENDENCIES |
|---|
| 50 | ! ------------ |
|---|
| 51 | use geometry, only: ngrid, nslope |
|---|
| 52 | |
|---|
| 53 | ! DECLARATION |
|---|
| 54 | ! ----------- |
|---|
| 55 | implicit none |
|---|
| 56 | |
|---|
| 57 | ! CODE |
|---|
| 58 | ! ---- |
|---|
| 59 | if (.not. allocated(albedo_PCM)) allocate(albedo_PCM(ngrid,nslope)) |
|---|
| 60 | if (.not. allocated(albedodat_PCM)) allocate(albedodat_PCM(ngrid)) |
|---|
| 61 | if (.not. allocated(emissivity_PCM)) allocate(emissivity_PCM(ngrid,nslope)) |
|---|
| 62 | |
|---|
| 63 | END SUBROUTINE ini_surface |
|---|
| 64 | !======================================================================= |
|---|
| 65 | |
|---|
| 66 | !======================================================================= |
|---|
| 67 | SUBROUTINE end_surface() |
|---|
| 68 | !----------------------------------------------------------------------- |
|---|
| 69 | ! NAME |
|---|
| 70 | ! end_surface |
|---|
| 71 | ! |
|---|
| 72 | ! DESCRIPTION |
|---|
| 73 | ! Deallocate surface arrays. |
|---|
| 74 | ! |
|---|
| 75 | ! AUTHORS & DATE |
|---|
| 76 | ! JB Clement, 12/2025 |
|---|
| 77 | ! |
|---|
| 78 | ! NOTES |
|---|
| 79 | ! |
|---|
| 80 | !----------------------------------------------------------------------- |
|---|
| 81 | |
|---|
| 82 | ! DECLARATION |
|---|
| 83 | ! ----------- |
|---|
| 84 | implicit none |
|---|
| 85 | |
|---|
| 86 | ! CODE |
|---|
| 87 | ! ---- |
|---|
| 88 | if (allocated(albedo_PCM)) deallocate(albedo_PCM) |
|---|
| 89 | if (allocated(albedodat_PCM)) deallocate(albedodat_PCM) |
|---|
| 90 | if (allocated(emissivity_PCM)) deallocate(emissivity_PCM ) |
|---|
| 91 | |
|---|
| 92 | END SUBROUTINE end_surface |
|---|
| 93 | !======================================================================= |
|---|
| 94 | |
|---|
| 95 | !======================================================================= |
|---|
| 96 | SUBROUTINE set_albedodat_PCM(albedodat_PCM_in) |
|---|
| 97 | !----------------------------------------------------------------------- |
|---|
| 98 | ! NAME |
|---|
| 99 | ! set_albedodat_PCM |
|---|
| 100 | ! |
|---|
| 101 | ! DESCRIPTION |
|---|
| 102 | ! Setter for 'albedodat_PCM'. |
|---|
| 103 | ! |
|---|
| 104 | ! AUTHORS & DATE |
|---|
| 105 | ! JB Clement, 12/2025 |
|---|
| 106 | ! |
|---|
| 107 | ! NOTES |
|---|
| 108 | ! |
|---|
| 109 | !----------------------------------------------------------------------- |
|---|
| 110 | |
|---|
| 111 | ! DECLARATION |
|---|
| 112 | ! ----------- |
|---|
| 113 | implicit none |
|---|
| 114 | |
|---|
| 115 | ! ARGUMENTS |
|---|
| 116 | ! --------- |
|---|
| 117 | real(dp), dimension(:), intent(in) :: albedodat_PCM_in |
|---|
| 118 | |
|---|
| 119 | ! CODE |
|---|
| 120 | ! ---- |
|---|
| 121 | albedodat_PCM(:) = albedodat_PCM_in(:) |
|---|
| 122 | |
|---|
| 123 | END SUBROUTINE set_albedodat_PCM |
|---|
| 124 | !======================================================================= |
|---|
| 125 | |
|---|
| 126 | !======================================================================= |
|---|
| 127 | SUBROUTINE set_albedo_PCM(albedo_PCM_in) |
|---|
| 128 | !----------------------------------------------------------------------- |
|---|
| 129 | ! NAME |
|---|
| 130 | ! set_albedo_PCM |
|---|
| 131 | ! |
|---|
| 132 | ! DESCRIPTION |
|---|
| 133 | ! Setter for 'albedo_PCM'. |
|---|
| 134 | ! |
|---|
| 135 | ! AUTHORS & DATE |
|---|
| 136 | ! JB Clement, 12/2025 |
|---|
| 137 | ! |
|---|
| 138 | ! NOTES |
|---|
| 139 | ! |
|---|
| 140 | !----------------------------------------------------------------------- |
|---|
| 141 | |
|---|
| 142 | ! DECLARATION |
|---|
| 143 | ! ----------- |
|---|
| 144 | implicit none |
|---|
| 145 | |
|---|
| 146 | ! ARGUMENTS |
|---|
| 147 | ! --------- |
|---|
| 148 | real(dp), dimension(:,:), intent(in) :: albedo_PCM_in |
|---|
| 149 | |
|---|
| 150 | ! CODE |
|---|
| 151 | ! ---- |
|---|
| 152 | albedo_PCM(:,:) = albedo_PCM_in(:,:) |
|---|
| 153 | |
|---|
| 154 | END SUBROUTINE set_albedo_PCM |
|---|
| 155 | !======================================================================= |
|---|
| 156 | |
|---|
| 157 | !======================================================================= |
|---|
| 158 | SUBROUTINE set_emissivity_PCM(emissivity_PCM_in) |
|---|
| 159 | !----------------------------------------------------------------------- |
|---|
| 160 | ! NAME |
|---|
| 161 | ! set_emissivity_PCM |
|---|
| 162 | ! |
|---|
| 163 | ! DESCRIPTION |
|---|
| 164 | ! Setter for 'emissivity_PCM'. |
|---|
| 165 | ! |
|---|
| 166 | ! AUTHORS & DATE |
|---|
| 167 | ! JB Clement, 12/2025 |
|---|
| 168 | ! |
|---|
| 169 | ! NOTES |
|---|
| 170 | ! |
|---|
| 171 | !----------------------------------------------------------------------- |
|---|
| 172 | |
|---|
| 173 | ! DECLARATION |
|---|
| 174 | ! ----------- |
|---|
| 175 | implicit none |
|---|
| 176 | |
|---|
| 177 | ! ARGUMENTS |
|---|
| 178 | ! --------- |
|---|
| 179 | real(dp), dimension(:,:), intent(in) :: emissivity_PCM_in |
|---|
| 180 | |
|---|
| 181 | ! CODE |
|---|
| 182 | ! ---- |
|---|
| 183 | emissivity_PCM(:,:) = emissivity_PCM_in(:,:) |
|---|
| 184 | |
|---|
| 185 | END SUBROUTINE set_emissivity_PCM |
|---|
| 186 | !======================================================================= |
|---|
| 187 | |
|---|
| 188 | !======================================================================= |
|---|
| 189 | SUBROUTINE build4PCM_surf_rad_prop(h2o_ice,co2_ice,albedo4PCM,emissivity4PCM) |
|---|
| 190 | !----------------------------------------------------------------------- |
|---|
| 191 | ! NAME |
|---|
| 192 | ! build4PCM_surf_rad_prop |
|---|
| 193 | ! |
|---|
| 194 | ! DESCRIPTION |
|---|
| 195 | ! Reconstructs albedo and emissivity for the PCM. |
|---|
| 196 | ! |
|---|
| 197 | ! AUTHORS & DATE |
|---|
| 198 | ! JB Clement, 12/2025 |
|---|
| 199 | ! C. Metz, 01/2026 |
|---|
| 200 | ! |
|---|
| 201 | ! NOTES |
|---|
| 202 | ! |
|---|
| 203 | !----------------------------------------------------------------------- |
|---|
| 204 | |
|---|
| 205 | ! DEPENDENCIES |
|---|
| 206 | ! ------------ |
|---|
| 207 | use geometry, only: ngrid, nslope, latitudes |
|---|
| 208 | use frost, only: h2o_frost4PCM, co2_frost4PCM |
|---|
| 209 | use display, only: print_msg |
|---|
| 210 | |
|---|
| 211 | ! DECLARATION |
|---|
| 212 | ! ----------- |
|---|
| 213 | implicit none |
|---|
| 214 | |
|---|
| 215 | ! ARGUMENTS |
|---|
| 216 | ! --------- |
|---|
| 217 | real(dp), dimension(:,:), intent(in) :: h2o_ice, co2_ice |
|---|
| 218 | real(dp), dimension(:,:), intent(out) :: albedo4PCM, emissivity4PCM |
|---|
| 219 | |
|---|
| 220 | ! LOCAL VARIABLES |
|---|
| 221 | ! --------------- |
|---|
| 222 | integer(di) :: i, islope, icap |
|---|
| 223 | real(dp) :: albedo_h2oice, albedo_h2ofrost, h2ofrost_albedo_threshold, emissivity_baresoil |
|---|
| 224 | real(dp), dimension(2) :: albedo_co2frost, emissivity_co2ice, albedo_co2ice |
|---|
| 225 | |
|---|
| 226 | ! CODE |
|---|
| 227 | ! ---- |
|---|
| 228 | ! Fetch parameters from "callphys.def" and "startfi.nc" |
|---|
| 229 | call read_albedo_emis(albedo_h2oice,albedo_h2ofrost,h2ofrost_albedo_threshold, & |
|---|
| 230 | albedo_co2ice,albedo_co2frost,emissivity_co2ice,emissivity_baresoil) |
|---|
| 231 | |
|---|
| 232 | ! Reconstruction Loop |
|---|
| 233 | call print_msg('> Building albedo and emmissivity for the PCM') |
|---|
| 234 | do i = 1,ngrid |
|---|
| 235 | ! Determine hemisphere: 1 = Northern, 2 = Southern |
|---|
| 236 | if (latitudes(i) < 0._dp) then |
|---|
| 237 | icap = 2 |
|---|
| 238 | else |
|---|
| 239 | icap = 1 |
|---|
| 240 | end if |
|---|
| 241 | |
|---|
| 242 | do islope = 1,nslope |
|---|
| 243 | ! Bare ground (default) |
|---|
| 244 | albedo4PCM(i,islope) = albedodat_PCM(i) |
|---|
| 245 | emissivity4PCM(i,islope) = emissivity_baresoil |
|---|
| 246 | |
|---|
| 247 | ! H2O ice |
|---|
| 248 | if (h2o_ice(i,islope) > 0._dp) then |
|---|
| 249 | albedo4PCM(i,islope) = albedo_h2oice |
|---|
| 250 | emissivity4PCM(i,islope) = 1._dp |
|---|
| 251 | end if |
|---|
| 252 | |
|---|
| 253 | ! CO2 ice (dominant over H2O ice) |
|---|
| 254 | if (co2_ice(i,islope) > 0._dp) then |
|---|
| 255 | albedo4PCM(i,islope) = albedo_co2ice(icap) |
|---|
| 256 | emissivity4PCM(i,islope) = albedo_co2ice(icap) |
|---|
| 257 | end if |
|---|
| 258 | |
|---|
| 259 | ! H2O frost (subject to threshold) |
|---|
| 260 | if (h2o_frost4PCM(i,islope) > h2ofrost_albedo_threshold) then |
|---|
| 261 | albedo4PCM(i,islope) = albedo_h2ofrost |
|---|
| 262 | emissivity4PCM(i,islope) = 1._dp |
|---|
| 263 | end if |
|---|
| 264 | |
|---|
| 265 | ! CO2 frost (final priority) |
|---|
| 266 | if (co2_frost4PCM(i,islope) > 0.) then |
|---|
| 267 | albedo4PCM(i,islope) = albedo_co2frost(icap) |
|---|
| 268 | emissivity4PCM(i,islope) = emissivity_co2ice(icap) |
|---|
| 269 | end if |
|---|
| 270 | end do |
|---|
| 271 | end do |
|---|
| 272 | |
|---|
| 273 | END SUBROUTINE build4PCM_surf_rad_prop |
|---|
| 274 | !======================================================================= |
|---|
| 275 | |
|---|
| 276 | !======================================================================= |
|---|
| 277 | SUBROUTINE read_albedo_emis(albedo_h2oice,albedo_h2ofrost,h2ofrost_albedo_threshold, & |
|---|
| 278 | albedo_co2ice,albedo_co2frost,emissivity_co2ice,emissivity_baresoil) |
|---|
| 279 | !----------------------------------------------------------------------- |
|---|
| 280 | ! NAME |
|---|
| 281 | ! read_albedo_emis |
|---|
| 282 | ! |
|---|
| 283 | ! DESCRIPTION |
|---|
| 284 | ! Reads albedo/emissivity parameters from "callphys.def" and |
|---|
| 285 | ! "startfi.nc" ('controle'). |
|---|
| 286 | ! |
|---|
| 287 | ! AUTHORS & DATE |
|---|
| 288 | ! C. Metz, 01/2026 |
|---|
| 289 | ! |
|---|
| 290 | ! NOTES |
|---|
| 291 | ! |
|---|
| 292 | !----------------------------------------------------------------------- |
|---|
| 293 | |
|---|
| 294 | ! DEPENDENCIES |
|---|
| 295 | ! ------------ |
|---|
| 296 | #ifdef CPP_IOIPSL |
|---|
| 297 | use IOIPSL, only: getin |
|---|
| 298 | #else |
|---|
| 299 | use ioipsl_getincom, only: getin |
|---|
| 300 | #endif |
|---|
| 301 | use config, only: callphys_name |
|---|
| 302 | use io_netcdf, only: open_nc, close_nc, startfi_name, get_dim_nc, get_var_nc |
|---|
| 303 | use stoppage, only: stop_clean |
|---|
| 304 | |
|---|
| 305 | ! DECLARATION |
|---|
| 306 | ! ----------- |
|---|
| 307 | implicit none |
|---|
| 308 | |
|---|
| 309 | ! ARGUMENTS |
|---|
| 310 | ! --------- |
|---|
| 311 | real(dp), intent(out) :: albedo_h2oice, albedo_h2ofrost, h2ofrost_albedo_threshold, emissivity_baresoil |
|---|
| 312 | real(dp), dimension(2), intent(out) :: albedo_co2frost, emissivity_co2ice, albedo_co2ice |
|---|
| 313 | |
|---|
| 314 | ! LOCAL VARIABLES |
|---|
| 315 | ! --------------- |
|---|
| 316 | logical(k4) :: here, cst_h2oice_albedo |
|---|
| 317 | integer(di) :: i, nindex |
|---|
| 318 | real(dp), dimension(:), allocatable :: controle |
|---|
| 319 | |
|---|
| 320 | ! CODE |
|---|
| 321 | ! ---- |
|---|
| 322 | inquire(file = callphys_name,exist = here) |
|---|
| 323 | if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find "'//callphys_name//'"!', 1) |
|---|
| 324 | |
|---|
| 325 | ! Read albedos of H2O frost and perennial ice from "callphys.def" |
|---|
| 326 | albedo_h2oice = 0.35_dp ! Default |
|---|
| 327 | call getin('albedo_h2o_cap',albedo_h2oice) |
|---|
| 328 | if (albedo_h2oice < 0._dp .or. albedo_h2oice > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_h2o_cap'' in "'//callphys_name//'" is out of bounds [0,1]!',1) |
|---|
| 329 | |
|---|
| 330 | albedo_h2ofrost = 0.35_dp ! Default |
|---|
| 331 | call getin('albedo_h2ofrost',albedo_h2ofrost) |
|---|
| 332 | if (albedo_h2ofrost < 0._dp .or. albedo_h2ofrost > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_h2ofrost'' in "'//callphys_name//'" is out of bounds [0,1]!',1) |
|---|
| 333 | |
|---|
| 334 | h2ofrost_albedo_threshold = 0.005_dp ! Default [kg/m2] |
|---|
| 335 | call getin('frost_albedo_threshold',h2ofrost_albedo_threshold) |
|---|
| 336 | |
|---|
| 337 | cst_h2oice_albedo = .false. ! Default |
|---|
| 338 | call getin('cst_cap_albedo',cst_h2oice_albedo) |
|---|
| 339 | if (cst_h2oice_albedo) albedo_h2ofrost = albedo_h2oice ! If true, then we don't account for water frost albedo effect |
|---|
| 340 | |
|---|
| 341 | ! Read albedos of CO2 perennial ice from "callphys.def" |
|---|
| 342 | albedo_co2ice(1) = 0.6_dp ! Default, north hemisphere |
|---|
| 343 | albedo_co2ice(2) = 0.85_dp ! Default, south hemisphere |
|---|
| 344 | call getin('albedo_co2ice_north',albedo_co2ice(1)) |
|---|
| 345 | call getin('albedo_co2ice_south',albedo_co2ice(2)) |
|---|
| 346 | if (albedo_co2ice(1) < 0._dp .or. albedo_co2ice(1) > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_co2ice_north'' in "'//callphys_name//'" is out of bounds [0,1]!',1) |
|---|
| 347 | if (albedo_co2ice(2) < 0._dp .or. albedo_co2ice(2) > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_co2ice_south'' in "'//callphys_name//'" is out of bounds [0,1]!',1) |
|---|
| 348 | |
|---|
| 349 | ! Read albedos of CO2 frost from "callphys.def" |
|---|
| 350 | call open_nc(startfi_name,'read') |
|---|
| 351 | call get_dim_nc('index',nindex) |
|---|
| 352 | allocate(controle(nindex)) |
|---|
| 353 | call get_var_nc('controle',controle) |
|---|
| 354 | call close_nc(startfi_name) |
|---|
| 355 | albedo_co2frost(1) = controle(22) |
|---|
| 356 | albedo_co2frost(2) = controle(23) |
|---|
| 357 | emissivity_co2ice(1) = controle(24) |
|---|
| 358 | emissivity_co2ice(2) = controle(25) |
|---|
| 359 | emissivity_baresoil = controle(26) |
|---|
| 360 | deallocate(controle) |
|---|
| 361 | do i = 1,2 |
|---|
| 362 | if (albedo_co2frost(i) < 0._dp .or. albedo_co2frost(i) > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_co2frost'' from ''controle(22:23)'' in "'//startfi_name//'" is out of bounds [0,1]!',1) |
|---|
| 363 | if (emissivity_co2ice(i) < 0._dp .or. emissivity_co2ice(i) > 1._dp) call stop_clean(__FILE__,__LINE__,'''emissivity_co2ice'' from ''controle(24:25)'' in "'//startfi_name//'" is out of bounds [0,1]!',1) |
|---|
| 364 | end do |
|---|
| 365 | if (emissivity_baresoil < 0._dp .or. emissivity_baresoil > 1._dp) call stop_clean(__FILE__,__LINE__,'''emissivity_baresoil'' from ''controle(26)'' in "'//startfi_name//'" is out of bounds [0,1]!',1) |
|---|
| 366 | |
|---|
| 367 | END SUBROUTINE read_albedo_emis |
|---|
| 368 | !======================================================================= |
|---|
| 369 | |
|---|
| 370 | END MODULE surface |
|---|