Line data Source code
1 : !> \file mo_edk_types.f90
2 : !> \brief \copybrief mo_edk_types
3 : !> \details \copydetails mo_edk_types
4 :
5 : !> \brief Module for EDK distance Type.
6 : !> \version 0.1
7 : !> \authors Sebastian Mueller
8 : !> \date 21.10.2020
9 : !> \details This file contains memory efficient types for EDK.
10 : !> \copyright Copyright 2005-\today, the CHS Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
11 : !! EDK is released under the LGPLv3+ license \license_note
12 : module mo_edk_types
13 : use mo_kind, only: i4, dp
14 : implicit none
15 : !> \class dist_t
16 : !> \brief Cell distances on demand.
17 : type dist_t
18 : integer(i4) :: ncell !< no of cells
19 : integer(i4) :: nstat !< no of stations
20 : real(dp) :: noDataValue !< masking value
21 : real(dp), dimension(:,:), allocatable :: cell_pos !< (x,y) pos of cells (i, j) - i: cell index; j: 1 - x, 2 - y
22 : real(dp), dimension(:,:), allocatable :: stat_pos !< (x,y) pos of stations (i, j) - i: station index; j: 1 - x, 2 - y
23 : real(dp), dimension(:,:), allocatable :: stat_z !< z values (i,j) - i: station; j - timeID
24 : contains
25 : !> \copydoc mo_edk_types::clean_dist
26 : procedure :: clean => clean_dist !< \see mo_edk_types::clean_dist
27 : !> \copydoc mo_edk_types::get_dist_cell_station
28 : procedure :: getCS => get_dist_cell_station !< \see mo_edk_types::get_dist_cell_station
29 : !> \copydoc mo_edk_types::get_dist_station_station
30 : procedure :: getSS => get_dist_station_station !< \see mo_edk_types::get_dist_station_station
31 : !> \copydoc mo_edk_types::get_squared_z_diff
32 : procedure :: getZ2 => get_squared_z_diff !< \see mo_edk_types::get_squared_z_diff
33 : end type dist_t
34 :
35 : contains
36 :
37 : !> \brief Clean up.
38 2 : subroutine clean_dist(self)
39 : implicit none
40 : class(dist_t), intent(inout) :: self
41 2 : if ( allocated(self%cell_pos)) deallocate (self%cell_pos)
42 2 : if ( allocated(self%stat_pos)) deallocate (self%stat_pos)
43 2 : if ( allocated(self%stat_z)) deallocate (self%stat_z)
44 2 : end subroutine clean_dist
45 :
46 : !> \brief Get distance between cell i and station j.
47 : !> \return Distance.
48 70248 : real(dp) function get_dist_cell_station(self, i, j) result(distCS)
49 : implicit none
50 : class(dist_t), intent(in) :: self
51 : integer, intent(in) :: i, j
52 70248 : distCS = dsqrt( ( self%cell_pos(i,1) - self%stat_pos(j,1) )**2 + ( self%cell_pos(i,2) - self%stat_pos(j,2) )** 2)
53 2 : end function get_dist_cell_station
54 :
55 : !> \brief Get distance between station i and station j.
56 : !> \return Distance.
57 642186 : real(dp) function get_dist_station_station(self, i, j) result(distSS)
58 : implicit none
59 : class(dist_t), intent(in) :: self
60 : integer, intent(in) :: i, j
61 642186 : distSS = dsqrt( ( self%stat_pos(i,1) - self%stat_pos(j,1) )**2 + ( self%stat_pos(i,2) - self%stat_pos(j,2) )** 2)
62 70248 : end function get_dist_station_station
63 :
64 : !> \brief Get squared elevation difference between station i and station j.
65 : !> \return Elevation difference squared.
66 106215 : real(dp) function get_squared_z_diff(self, i, j, jd) result(distZ2)
67 : implicit none
68 : class(dist_t), intent(in) :: self
69 : integer, intent(in) :: i, j, jd
70 212430 : if ( (self%stat_z(i, jd) /= self%noDataValue ) .and. &
71 106215 : (self%stat_z(j, jd) /= self%noDataValue ) ) then
72 87600 : distZ2 = ( self%stat_z(i, jd) - self%stat_z(j, jd) )**2
73 : else
74 : distZ2 = self%noDataValue
75 : end if
76 642186 : end function get_squared_z_diff
77 :
78 106215 : end module mo_edk_types
|