| 1 | MODULE geometry |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! geometry |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! Grid metrics/parameters and basic tools to convert between |
|---|
| 8 | ! different grids. |
|---|
| 9 | ! |
|---|
| 10 | ! AUTHORS & DATE |
|---|
| 11 | ! JB Clement, 12/2025 |
|---|
| 12 | ! |
|---|
| 13 | ! NOTES |
|---|
| 14 | ! |
|---|
| 15 | !----------------------------------------------------------------------- |
|---|
| 16 | |
|---|
| 17 | ! DEPENDENCIES |
|---|
| 18 | ! ------------ |
|---|
| 19 | use numerics, only: dp, di, k4 |
|---|
| 20 | |
|---|
| 21 | ! DECLARATION |
|---|
| 22 | ! ----------- |
|---|
| 23 | implicit none |
|---|
| 24 | |
|---|
| 25 | ! PARAMETERS |
|---|
| 26 | ! ---------- |
|---|
| 27 | integer(di), protected :: nlon ! Number of longitudes |
|---|
| 28 | integer(di), protected :: nlat ! Number of latitudes |
|---|
| 29 | integer(di), protected :: ngrid ! Number of grid points |
|---|
| 30 | integer(di), protected :: nlayer ! Number of atmospheric layers |
|---|
| 31 | integer(di), protected :: nsoil_PCM ! Number of soil layers in the PCM |
|---|
| 32 | integer(di), parameter :: nsoil = 68_di ! Number of soil layers in the PEM |
|---|
| 33 | integer(di), protected :: nslope ! Number of slopes |
|---|
| 34 | integer(di), protected :: nday ! Number of sols in a year |
|---|
| 35 | real(dp), dimension(:), allocatable, protected :: longitudes ! Longitudes |
|---|
| 36 | real(dp), dimension(:), allocatable, protected :: latitudes ! Latitudes |
|---|
| 37 | real(dp), dimension(:), allocatable, protected :: cell_area ! Cell area |
|---|
| 38 | real(dp), protected :: total_surface ! Total surface of the grid/planet |
|---|
| 39 | logical(k4), protected :: dim_init = .false. ! Flag to true if dimensions are well initialized |
|---|
| 40 | |
|---|
| 41 | contains |
|---|
| 42 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 43 | |
|---|
| 44 | !======================================================================= |
|---|
| 45 | SUBROUTINE ini_geometry() |
|---|
| 46 | !----------------------------------------------------------------------- |
|---|
| 47 | ! NAME |
|---|
| 48 | ! ini_geometry |
|---|
| 49 | ! |
|---|
| 50 | ! DESCRIPTION |
|---|
| 51 | ! Initialize the parameters of module 'geometry'. |
|---|
| 52 | ! |
|---|
| 53 | ! AUTHORS & DATE |
|---|
| 54 | ! JB Clement, 12/2025 |
|---|
| 55 | ! |
|---|
| 56 | ! NOTES |
|---|
| 57 | ! |
|---|
| 58 | !----------------------------------------------------------------------- |
|---|
| 59 | |
|---|
| 60 | ! DEPENDENCIES |
|---|
| 61 | ! ------------ |
|---|
| 62 | use stoppage, only: stop_clean |
|---|
| 63 | use io_netcdf, only: startfi_name, xios_day_name1, open_nc, close_nc, get_dim_nc, get_var_nc |
|---|
| 64 | use display, only: print_msg |
|---|
| 65 | use utility, only: int2str |
|---|
| 66 | |
|---|
| 67 | ! DECLARATION |
|---|
| 68 | ! ----------- |
|---|
| 69 | implicit none |
|---|
| 70 | |
|---|
| 71 | ! LOCAL VARIABLES |
|---|
| 72 | ! --------------- |
|---|
| 73 | logical(k4) :: here, found |
|---|
| 74 | |
|---|
| 75 | ! CODE |
|---|
| 76 | ! ---- |
|---|
| 77 | ! Reading from "startfi.nc" file |
|---|
| 78 | inquire(file = startfi_name,exist = here) |
|---|
| 79 | if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//startfi_name//'"!',1) |
|---|
| 80 | call open_nc(startfi_name,'read') |
|---|
| 81 | call get_dim_nc('physical_points',ngrid) |
|---|
| 82 | call get_dim_nc('subsurface_layers',nsoil_PCM) |
|---|
| 83 | call get_dim_nc('nlayer',nlayer) |
|---|
| 84 | call get_dim_nc('nslope',nslope,found) |
|---|
| 85 | if (.not. found) nslope = 1 |
|---|
| 86 | if (.not. allocated(longitudes)) allocate(longitudes(ngrid)) |
|---|
| 87 | if (.not. allocated(latitudes)) allocate(latitudes(ngrid)) |
|---|
| 88 | if (.not. allocated(cell_area)) allocate(cell_area(ngrid)) |
|---|
| 89 | call get_var_nc('longitude',longitudes) |
|---|
| 90 | call get_var_nc('latitude',latitudes) |
|---|
| 91 | call get_var_nc('area',cell_area) |
|---|
| 92 | call close_nc(startfi_name) |
|---|
| 93 | |
|---|
| 94 | ! Reading from XIOS daily output file |
|---|
| 95 | inquire(file = xios_day_name1,exist = here) |
|---|
| 96 | if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//xios_day_name1//'"!',1) |
|---|
| 97 | call open_nc(xios_day_name1,'read') |
|---|
| 98 | call get_dim_nc('lon',nlon) |
|---|
| 99 | call get_dim_nc('lat',nlat) |
|---|
| 100 | call get_dim_nc('time_counter',nday) |
|---|
| 101 | call close_nc(xios_day_name1) |
|---|
| 102 | |
|---|
| 103 | ! Total surface |
|---|
| 104 | total_surface = sum(cell_area) |
|---|
| 105 | |
|---|
| 106 | ! State that dimentions are well initialized |
|---|
| 107 | dim_init = .true. |
|---|
| 108 | call print_msg('> Reading dimensions') |
|---|
| 109 | call print_msg('nlon = '//int2str(nlon)) |
|---|
| 110 | call print_msg('nlat = '//int2str(nlat)) |
|---|
| 111 | call print_msg('ngrid = '//int2str(ngrid)) |
|---|
| 112 | call print_msg('nslope = '//int2str(nslope)) |
|---|
| 113 | call print_msg('nlayer = '//int2str(nlayer)) |
|---|
| 114 | call print_msg('nsoil_PCM = '//int2str(nsoil_PCM)) |
|---|
| 115 | call print_msg('nsoil = '//int2str(nsoil)) |
|---|
| 116 | call print_msg('nday = '//int2str(nday)) |
|---|
| 117 | |
|---|
| 118 | ! Check consistency of grids |
|---|
| 119 | if (ngrid /= 2 + nlon*(nlat - 2)) call stop_clean(__FILE__,__LINE__,'geometry inconsistency between nlon/nlat and ngrid!',1) |
|---|
| 120 | if (nsoil_PCM > nsoil) call stop_clean(__FILE__,__LINE__,'geometry inconsistency between nsoil and nsoil_PCM! nsoil must be greater than nsoil_PCM.',1) |
|---|
| 121 | |
|---|
| 122 | END SUBROUTINE ini_geometry |
|---|
| 123 | !======================================================================= |
|---|
| 124 | |
|---|
| 125 | !======================================================================= |
|---|
| 126 | SUBROUTINE end_geometry() |
|---|
| 127 | !----------------------------------------------------------------------- |
|---|
| 128 | ! NAME |
|---|
| 129 | ! end_geometry |
|---|
| 130 | ! |
|---|
| 131 | ! DESCRIPTION |
|---|
| 132 | ! Deallocate geometry arrays. |
|---|
| 133 | ! |
|---|
| 134 | ! AUTHORS & DATE |
|---|
| 135 | ! JB Clement, 12/2025 |
|---|
| 136 | ! |
|---|
| 137 | ! NOTES |
|---|
| 138 | ! |
|---|
| 139 | !----------------------------------------------------------------------- |
|---|
| 140 | |
|---|
| 141 | ! DECLARATION |
|---|
| 142 | ! ----------- |
|---|
| 143 | implicit none |
|---|
| 144 | |
|---|
| 145 | ! CODE |
|---|
| 146 | ! ---- |
|---|
| 147 | if (allocated(longitudes)) deallocate(longitudes) |
|---|
| 148 | if (allocated(latitudes)) deallocate(latitudes) |
|---|
| 149 | if (allocated(cell_area)) deallocate(cell_area) |
|---|
| 150 | |
|---|
| 151 | END SUBROUTINE end_geometry |
|---|
| 152 | !======================================================================= |
|---|
| 153 | |
|---|
| 154 | !======================================================================= |
|---|
| 155 | SUBROUTINE lonlat2vect(nlon,nlat,ngrid,v_ll,v_vect) |
|---|
| 156 | !----------------------------------------------------------------------- |
|---|
| 157 | ! NAME |
|---|
| 158 | ! lonlat2vect |
|---|
| 159 | ! |
|---|
| 160 | ! DESCRIPTION |
|---|
| 161 | ! Convert data from lon x lat grid (where values at the poles are |
|---|
| 162 | ! duplicated) into a vector. |
|---|
| 163 | ! |
|---|
| 164 | ! AUTHORS & DATE |
|---|
| 165 | ! JB Clement, 12/2025 |
|---|
| 166 | ! |
|---|
| 167 | ! NOTES |
|---|
| 168 | ! The longitudes -180 and +180 are not duplicated like in the PCM |
|---|
| 169 | ! dynamics. |
|---|
| 170 | !----------------------------------------------------------------------- |
|---|
| 171 | |
|---|
| 172 | ! DEPENDENCIES |
|---|
| 173 | ! ------------ |
|---|
| 174 | use stoppage, only: stop_clean |
|---|
| 175 | |
|---|
| 176 | ! DECLARATION |
|---|
| 177 | ! ----------- |
|---|
| 178 | implicit none |
|---|
| 179 | |
|---|
| 180 | ! ARGUMENTS |
|---|
| 181 | !---------- |
|---|
| 182 | integer(di), intent(in) :: nlon, nlat, ngrid |
|---|
| 183 | real(dp), dimension(nlon,nlat), intent(in) :: v_ll |
|---|
| 184 | real(dp), dimension(ngrid), intent(out) :: v_vect |
|---|
| 185 | |
|---|
| 186 | ! LOCAL VARIABLES |
|---|
| 187 | !---------------- |
|---|
| 188 | integer(di) :: i, j |
|---|
| 189 | |
|---|
| 190 | ! CODE |
|---|
| 191 | !----- |
|---|
| 192 | if (ngrid == 1) then ! 1D case |
|---|
| 193 | v_vect(1) = v_ll(1,1) |
|---|
| 194 | return |
|---|
| 195 | else ! 3D |
|---|
| 196 | ! Check |
|---|
| 197 | if (ngrid /= nlon*(nlat - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1) |
|---|
| 198 | |
|---|
| 199 | ! Initialization |
|---|
| 200 | v_vect = 0. |
|---|
| 201 | |
|---|
| 202 | ! Treatment of the poles |
|---|
| 203 | v_vect(1) = v_ll(1,1) |
|---|
| 204 | v_vect(ngrid) = v_ll(1,nlat) |
|---|
| 205 | |
|---|
| 206 | ! Treatment of regular points |
|---|
| 207 | do j = 2,nlat - 1 |
|---|
| 208 | do i = 1,nlon |
|---|
| 209 | v_vect(1 + i + (j - 2)*nlon) = v_ll(i,j) |
|---|
| 210 | end do |
|---|
| 211 | end do |
|---|
| 212 | end if |
|---|
| 213 | |
|---|
| 214 | END SUBROUTINE lonlat2vect |
|---|
| 215 | !======================================================================= |
|---|
| 216 | |
|---|
| 217 | !======================================================================= |
|---|
| 218 | SUBROUTINE vect2lonlat(nlon,nlat,ngrid,v_vect,v_ll) |
|---|
| 219 | !----------------------------------------------------------------------- |
|---|
| 220 | ! NAME |
|---|
| 221 | ! vect2lonlat |
|---|
| 222 | ! |
|---|
| 223 | ! DESCRIPTION |
|---|
| 224 | ! Convert data from a vector into lon x lat grid (where values |
|---|
| 225 | ! at the poles are duplicated). |
|---|
| 226 | ! |
|---|
| 227 | ! AUTHORS & DATE |
|---|
| 228 | ! JB Clement, 12/2025 |
|---|
| 229 | ! |
|---|
| 230 | ! NOTES |
|---|
| 231 | ! The longitudes -180 and +180 are not duplicated like in the PCM |
|---|
| 232 | ! dynamics. |
|---|
| 233 | !----------------------------------------------------------------------- |
|---|
| 234 | |
|---|
| 235 | ! DEPENDENCIES |
|---|
| 236 | ! ------------ |
|---|
| 237 | use stoppage, only: stop_clean |
|---|
| 238 | |
|---|
| 239 | ! DECLARATION |
|---|
| 240 | ! ----------- |
|---|
| 241 | implicit none |
|---|
| 242 | |
|---|
| 243 | ! ARGUMENTS |
|---|
| 244 | !---------- |
|---|
| 245 | integer(di), intent(in) :: nlon, nlat, ngrid |
|---|
| 246 | real(dp), dimension(ngrid), intent(in) :: v_vect |
|---|
| 247 | real(dp), dimension(nlon,nlat), intent(out) :: v_ll |
|---|
| 248 | |
|---|
| 249 | ! LOCAL VARIABLES |
|---|
| 250 | !---------------- |
|---|
| 251 | integer(di) :: i, j |
|---|
| 252 | |
|---|
| 253 | ! CODE |
|---|
| 254 | !----- |
|---|
| 255 | if (ngrid == 1) then ! 1D case |
|---|
| 256 | v_ll(1,1) = v_vect(1) |
|---|
| 257 | return |
|---|
| 258 | else ! 3D |
|---|
| 259 | ! Check |
|---|
| 260 | if (ngrid /= nlon*(nlat - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1) |
|---|
| 261 | |
|---|
| 262 | ! Initialization |
|---|
| 263 | v_ll = 0. |
|---|
| 264 | |
|---|
| 265 | ! Treatment of the poles |
|---|
| 266 | v_ll(:,1) = v_vect(1) |
|---|
| 267 | v_ll(:,nlat) = v_vect(ngrid) |
|---|
| 268 | |
|---|
| 269 | ! Treatment of regular points |
|---|
| 270 | do j = 2,nlat - 1 |
|---|
| 271 | do i = 1,nlon |
|---|
| 272 | v_ll(i,j) = v_vect(1 + i + (j - 2)*nlon) |
|---|
| 273 | end do |
|---|
| 274 | end do |
|---|
| 275 | end if |
|---|
| 276 | |
|---|
| 277 | END SUBROUTINE vect2lonlat |
|---|
| 278 | !======================================================================= |
|---|
| 279 | |
|---|
| 280 | !======================================================================= |
|---|
| 281 | SUBROUTINE dyngrd2vect(ni,nj,ngrid,v_dyn,v_vect) |
|---|
| 282 | !----------------------------------------------------------------------- |
|---|
| 283 | ! NAME |
|---|
| 284 | ! dyngrd2vect |
|---|
| 285 | ! |
|---|
| 286 | ! DESCRIPTION |
|---|
| 287 | ! Convert data from dynamical lon x lat grid (where values at the |
|---|
| 288 | ! poles are duplicated and longitude -180°/+180° is duplicated) into |
|---|
| 289 | ! a vector. |
|---|
| 290 | ! |
|---|
| 291 | ! AUTHORS & DATE |
|---|
| 292 | ! JB Clement, 12/2025 |
|---|
| 293 | ! |
|---|
| 294 | ! NOTES |
|---|
| 295 | ! The longitudes -180 and +180 are duplicated (PCM dynamics). |
|---|
| 296 | !----------------------------------------------------------------------- |
|---|
| 297 | |
|---|
| 298 | ! DEPENDENCIES |
|---|
| 299 | ! ------------ |
|---|
| 300 | use stoppage, only: stop_clean |
|---|
| 301 | |
|---|
| 302 | ! DECLARATION |
|---|
| 303 | ! ----------- |
|---|
| 304 | implicit none |
|---|
| 305 | |
|---|
| 306 | ! ARGUMENTS |
|---|
| 307 | !---------- |
|---|
| 308 | integer(di), intent(in) :: ni, nj, ngrid |
|---|
| 309 | real(dp), dimension(ni,nj), intent(in) :: v_dyn |
|---|
| 310 | real(dp), dimension(ngrid), intent(out) :: v_vect |
|---|
| 311 | |
|---|
| 312 | ! LOCAL VARIABLES |
|---|
| 313 | !---------------- |
|---|
| 314 | integer(di) :: i, j |
|---|
| 315 | |
|---|
| 316 | ! CODE |
|---|
| 317 | !----- |
|---|
| 318 | if (ngrid == 1) then ! 1D case |
|---|
| 319 | v_vect(1) = v_dyn(1,1) |
|---|
| 320 | return |
|---|
| 321 | else ! 3D |
|---|
| 322 | ! Check |
|---|
| 323 | if (ngrid /= (ni - 1)*(nj - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1) |
|---|
| 324 | |
|---|
| 325 | ! Initialization |
|---|
| 326 | v_vect = 0. |
|---|
| 327 | |
|---|
| 328 | ! Treatment of the poles |
|---|
| 329 | v_vect(1) = v_dyn(1,1) |
|---|
| 330 | v_vect(ngrid) = v_dyn(1,nj) |
|---|
| 331 | |
|---|
| 332 | ! Treatment of regular points |
|---|
| 333 | do j = 2,nj - 1 |
|---|
| 334 | do i = 1,ni - 1 |
|---|
| 335 | v_vect(1 + i + (j - 2)*(ni - 1)) = v_dyn(i,j) |
|---|
| 336 | end do |
|---|
| 337 | end do |
|---|
| 338 | end if |
|---|
| 339 | |
|---|
| 340 | END SUBROUTINE dyngrd2vect |
|---|
| 341 | !======================================================================= |
|---|
| 342 | |
|---|
| 343 | !======================================================================= |
|---|
| 344 | SUBROUTINE vect2dyngrd(ni,nj,ngrid,v_vect,v_dyn) |
|---|
| 345 | !----------------------------------------------------------------------- |
|---|
| 346 | ! NAME |
|---|
| 347 | ! vect2dyngrd |
|---|
| 348 | ! |
|---|
| 349 | ! DESCRIPTION |
|---|
| 350 | ! Convert data from a vector into lon x lat grid (where values at the |
|---|
| 351 | ! poles are duplicated and longitude -180°/+180° is duplicated). |
|---|
| 352 | ! |
|---|
| 353 | ! AUTHORS & DATE |
|---|
| 354 | ! JB Clement, 12/2025 |
|---|
| 355 | ! |
|---|
| 356 | ! NOTES |
|---|
| 357 | ! The longitudes -180 and +180 are not duplicated like in the PCM |
|---|
| 358 | ! dynamics. |
|---|
| 359 | !----------------------------------------------------------------------- |
|---|
| 360 | |
|---|
| 361 | ! DEPENDENCIES |
|---|
| 362 | ! ------------ |
|---|
| 363 | use stoppage, only: stop_clean |
|---|
| 364 | |
|---|
| 365 | ! DECLARATION |
|---|
| 366 | ! ----------- |
|---|
| 367 | implicit none |
|---|
| 368 | |
|---|
| 369 | ! ARGUMENTS |
|---|
| 370 | !---------- |
|---|
| 371 | integer(di), intent(in) :: ni, nj, ngrid |
|---|
| 372 | real(dp), dimension(ngrid), intent(in) :: v_vect |
|---|
| 373 | real(dp), dimension(ni,nj), intent(out) :: v_dyn |
|---|
| 374 | |
|---|
| 375 | ! LOCAL VARIABLES |
|---|
| 376 | !---------------- |
|---|
| 377 | integer(di) :: i, j |
|---|
| 378 | |
|---|
| 379 | ! CODE |
|---|
| 380 | !----- |
|---|
| 381 | if (ngrid == 1) then ! 1D case |
|---|
| 382 | v_dyn(1,1) = v_vect(1) |
|---|
| 383 | return |
|---|
| 384 | else ! 3D |
|---|
| 385 | ! Check |
|---|
| 386 | if (ngrid /= (ni - 1)*(nj - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1) |
|---|
| 387 | |
|---|
| 388 | ! Initialization |
|---|
| 389 | v_dyn = 0. |
|---|
| 390 | |
|---|
| 391 | ! Treatment of the poles |
|---|
| 392 | v_dyn(:,1) = v_vect(1) |
|---|
| 393 | v_dyn(:,nj) = v_vect(ngrid) |
|---|
| 394 | |
|---|
| 395 | ! Treatment of regular points |
|---|
| 396 | do j = 2,nj - 1 |
|---|
| 397 | do i = 1,ni - 1 |
|---|
| 398 | v_dyn(i,j) = v_vect(1 + i + (j - 2)*(ni - 1)) |
|---|
| 399 | end do |
|---|
| 400 | v_dyn(ni,j) = v_dyn(1,j) |
|---|
| 401 | end do |
|---|
| 402 | end if |
|---|
| 403 | |
|---|
| 404 | END SUBROUTINE vect2dyngrd |
|---|
| 405 | !======================================================================= |
|---|
| 406 | |
|---|
| 407 | END MODULE geometry |
|---|