| 1 | MODULE surf_temp |
|---|
| 2 | |
|---|
| 3 | implicit none |
|---|
| 4 | |
|---|
| 5 | !======================================================================= |
|---|
| 6 | contains |
|---|
| 7 | !======================================================================= |
|---|
| 8 | |
|---|
| 9 | SUBROUTINE update_tsurf_nearest_baresoil(ngrid,nslope,nlon,nlat,latitude,tsurf_avg,co2_ice,is_co2ice_ini,co2ice_disappeared) |
|---|
| 10 | |
|---|
| 11 | use grid_conversion, only: vect2lonlat, lonlat2vect |
|---|
| 12 | |
|---|
| 13 | implicit none |
|---|
| 14 | |
|---|
| 15 | ! Inputs: |
|---|
| 16 | integer, intent(in) :: nlon, nlat, nslope, ngrid |
|---|
| 17 | real, dimension(ngrid,nslope), intent(in) :: co2_ice |
|---|
| 18 | real, dimension(ngrid), intent(in) :: latitude |
|---|
| 19 | logical, dimension(ngrid,nslope), intent(in) :: is_co2ice_ini |
|---|
| 20 | ! Outputs: |
|---|
| 21 | real, dimension(ngrid,nslope), intent(inout) :: tsurf_avg |
|---|
| 22 | logical, dimension(ngrid,nslope), intent(inout) :: co2ice_disappeared |
|---|
| 23 | ! Local variables: |
|---|
| 24 | real, parameter :: eps = 1.e-10 |
|---|
| 25 | integer :: islope, i, j, k, radius, rmax, di, dj, ii, jj |
|---|
| 26 | logical :: found |
|---|
| 27 | real, dimension(nlon,nlat,nslope) :: tsurf_ll, co2ice_ll, mask_co2ice_ini, co2ice_disappeared_ll |
|---|
| 28 | real, dimension(nlon,nlat) :: latitude_ll |
|---|
| 29 | real, dimension(ngrid) :: tmp |
|---|
| 30 | integer, dimension(nslope - 1) :: priority |
|---|
| 31 | |
|---|
| 32 | ! Check to escape the subroutine (not relevant in 1D) |
|---|
| 33 | if (ngrid == 1) return |
|---|
| 34 | |
|---|
| 35 | write(*,*) "> Updating surface temperature where ice disappeared" |
|---|
| 36 | ! Convert from reduced grid to lon-lat grid |
|---|
| 37 | call vect2lonlat(nlon,nlat,ngrid,latitude,latitude_ll) |
|---|
| 38 | do islope = 1,nslope |
|---|
| 39 | call vect2lonlat(nlon,nlat,ngrid,tsurf_avg(:,islope),tsurf_ll(:,:,islope)) |
|---|
| 40 | call vect2lonlat(nlon,nlat,ngrid,co2_ice(:,islope),co2ice_ll(:,:,islope)) |
|---|
| 41 | call vect2lonlat(nlon,nlat,ngrid,merge(1.,0.,is_co2ice_ini(:,islope)),mask_co2ice_ini(:,:,islope)) |
|---|
| 42 | call vect2lonlat(nlon,nlat,ngrid,merge(1.,0.,co2ice_disappeared(:,islope)),co2ice_disappeared_ll(:,:,islope)) |
|---|
| 43 | enddo |
|---|
| 44 | |
|---|
| 45 | ! For each point where ice disappeared |
|---|
| 46 | rmax = max(nlon,nlat) |
|---|
| 47 | do j = 1,nlat |
|---|
| 48 | do i = 1,nlon |
|---|
| 49 | do islope = 1,nslope |
|---|
| 50 | if (mask_co2ice_ini(i,j,islope) > 0.5 .and. co2ice_ll(i,j,islope) < eps .and. co2ice_disappeared_ll(i,j,islope) < 0.5) then |
|---|
| 51 | found = .false. |
|---|
| 52 | co2ice_disappeared_ll(i,j,islope) = 1. |
|---|
| 53 | call get_slope_priority(latitude_ll(i,j),nslope,islope,priority) |
|---|
| 54 | do k = 1,nslope - 1 |
|---|
| 55 | if (mask_co2ice_ini(i,j,priority(k)) < 0.5) then |
|---|
| 56 | tsurf_ll(i,j,islope) = tsurf_ll(i,j,priority(k)) |
|---|
| 57 | found = .true. |
|---|
| 58 | exit |
|---|
| 59 | endif |
|---|
| 60 | enddo |
|---|
| 61 | |
|---|
| 62 | radius = 1 |
|---|
| 63 | do while (.not. found .and. radius <= rmax) ! Only if no adjacent slopes holds bare soil |
|---|
| 64 | do dj = -radius,radius |
|---|
| 65 | do di = -radius,radius |
|---|
| 66 | if (abs(di) + abs(dj) == radius) then |
|---|
| 67 | ii = i + di |
|---|
| 68 | jj = j + dj |
|---|
| 69 | ! Longitudinal periodicity |
|---|
| 70 | if (ii < 1) then |
|---|
| 71 | ii = ii + nlon |
|---|
| 72 | else if (ii > nlon) then |
|---|
| 73 | ii = ii - nlon |
|---|
| 74 | endif |
|---|
| 75 | ! Latitude boundaries |
|---|
| 76 | if (jj >= 1 .and. jj <= nlat) then |
|---|
| 77 | call get_slope_priority(latitude_ll(ii,jj),nslope,islope,priority) |
|---|
| 78 | do k = 1,nslope - 1 |
|---|
| 79 | if (mask_co2ice_ini(ii,jj,priority(k)) < 0.5) then |
|---|
| 80 | tsurf_ll(i,j,islope) = tsurf_ll(ii,jj,priority(k)) |
|---|
| 81 | found = .true. |
|---|
| 82 | exit |
|---|
| 83 | endif |
|---|
| 84 | enddo |
|---|
| 85 | endif |
|---|
| 86 | endif |
|---|
| 87 | if (found) exit |
|---|
| 88 | enddo |
|---|
| 89 | if (found) exit |
|---|
| 90 | enddo |
|---|
| 91 | radius = radius + 1 |
|---|
| 92 | enddo |
|---|
| 93 | if (.not. found) write(*,*) "WARNING: no bare soil found for ice disappeared on point:",i,j,islope |
|---|
| 94 | endif |
|---|
| 95 | enddo |
|---|
| 96 | enddo |
|---|
| 97 | enddo |
|---|
| 98 | |
|---|
| 99 | ! Convert back from lon-lat grid to reduced grid |
|---|
| 100 | do islope = 1,nslope |
|---|
| 101 | call lonlat2vect(nlon,nlat,ngrid,tsurf_ll(:,:,islope),tsurf_avg(:,islope)) |
|---|
| 102 | call lonlat2vect(nlon,nlat,ngrid,co2ice_disappeared_ll(:,:,islope),tmp) |
|---|
| 103 | where (tmp > 0.5) co2ice_disappeared(:,islope) = .true. |
|---|
| 104 | enddo |
|---|
| 105 | |
|---|
| 106 | END SUBROUTINE update_tsurf_nearest_baresoil |
|---|
| 107 | |
|---|
| 108 | !======================================================================= |
|---|
| 109 | SUBROUTINE get_slope_priority(lat,nslope,islope,priority) |
|---|
| 110 | ! Priority given to equator-ward slope which are most likely to hold no ice |
|---|
| 111 | |
|---|
| 112 | implicit none |
|---|
| 113 | |
|---|
| 114 | ! Inputs: |
|---|
| 115 | real, intent(in) :: lat |
|---|
| 116 | integer, intent(in) :: nslope, islope |
|---|
| 117 | ! Outputs: |
|---|
| 118 | integer, dimension(nslope - 1), intent(out) :: priority |
|---|
| 119 | ! Locals: |
|---|
| 120 | integer :: i, k |
|---|
| 121 | |
|---|
| 122 | ! Code |
|---|
| 123 | !----- |
|---|
| 124 | k = 1 |
|---|
| 125 | |
|---|
| 126 | ! Northern hemisphere |
|---|
| 127 | if (lat > 0.) then |
|---|
| 128 | ! Equator-ward slopes |
|---|
| 129 | do i = islope - 1,1,-1 |
|---|
| 130 | priority(k) = i |
|---|
| 131 | k = k + 1 |
|---|
| 132 | enddo |
|---|
| 133 | ! Pole-ward slopes |
|---|
| 134 | do i = islope + 1,nslope |
|---|
| 135 | priority(k) = i |
|---|
| 136 | k = k + 1 |
|---|
| 137 | enddo |
|---|
| 138 | else ! Southern hemisphere |
|---|
| 139 | ! Equator-ward slopes |
|---|
| 140 | do i = islope + 1,nslope |
|---|
| 141 | priority(k) = i |
|---|
| 142 | k = k + 1 |
|---|
| 143 | enddo |
|---|
| 144 | ! Pole-ward slopes |
|---|
| 145 | do i = islope - 1,1,-1 |
|---|
| 146 | priority(k) = i |
|---|
| 147 | k = k + 1 |
|---|
| 148 | enddo |
|---|
| 149 | endif |
|---|
| 150 | |
|---|
| 151 | END SUBROUTINE get_slope_priority |
|---|
| 152 | |
|---|
| 153 | END MODULE surf_temp |
|---|