5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_read_latlon.f90
Go to the documentation of this file.
1!> \file mo_read_latlon.f90
2!> \brief \copybrief mo_read_latlon
3!> \details \copydetails mo_read_latlon
4
5!> \brief reading latitude and longitude coordinates for each domain
6!> \details This module provides routines for reading latitude and longitude coordinates from file.
7!> \authors Stephan Thober
8!> \date Nov 2013
9!> \copyright Copyright 2005-\today, the mHM Developers, Luis Samaniego, Sabine Attinger: All rights reserved.
10!! mHM is released under the LGPLv3+ license \license_note
11!> \ingroup f_common
13
14 USE mo_kind, ONLY : i4, dp
15 use mo_message, only: error_message
16 use mo_string_utils, only : num2str
17
18 ! Of course
19 IMPLICIT NONE
20
21 PUBLIC :: read_latlon
22
23 PRIVATE
24
25CONTAINS
26
27 ! ------------------------------------------------------------------
28
29 ! NAME
30 ! read_latlon
31
32 ! PURPOSE
33 !> \brief reads latitude and longitude coordinates
34
35 !> \details reads latitude and longitude coordinates from
36 !> netcdf file for each domain and appends it to the global
37 !> variables latitude and longitude.
38
39 ! INTENT(IN)
40 !> \param[in] "integer(i4) :: ii" domain indexFile name of the domains must be xxx_latlon.nc, wherexxx
41 !> is the domain id. Variable names in the netcdf filehave to be 'lat' for latitude and 'lon' for longitude.
42 !> \param[in] "character(*) :: lon_var_name"
43 !> \param[in] "character(*) :: lat_var_name"
44 !> \param[in] "character(*) :: level_name"
45
46 ! INTENT(INOUT)
47 !> \param[inout] "type(Grid) :: level"
48
49 ! HISTORY
50 !> \authors Stephan Thober
51
52 !> \date Nov 2013
53
54 ! Modifications:
55 ! Stephan Thober, Sep 2015 - added latitude and longitude for level 0
56 ! Stephan Thober, Oct 2015 - added L1_rect_latitude and L1_rect_longitude
57 ! David Schaefer, May 2016 - removed ncread dependency
58 ! Robert Schweppe, Mar 2018 - major rewrite
59 ! Robert Schweppe Jun 2018 - refactoring and reformatting
60
61 subroutine read_latlon(ii, lon_var_name, lat_var_name, level_name, level)
62
63 use mo_common_types, only: grid
65 use mo_netcdf, only : ncdataset, ncvariable
66
67 implicit none
68
69 ! domain indexFile name of the domains must be xxx_latlon.nc, wherexxx is the domain id. Variable names in the netcdf
70 ! filehave to be 'lat' for latitude and 'lon' for longitude.
71 integer(i4), intent(in) :: ii
72
73 character(*), intent(in) :: lon_var_name
74
75 character(*), intent(in) :: lat_var_name
76
77 character(*), intent(in) :: level_name
78
79 type(grid), intent(inout) :: level
80
81 ! file name
82 character(256) :: fname
83
84 ! dummy variable
85 real(dp), dimension(:, :), allocatable :: dummy
86
87 type(ncdataset) :: nc
88
89 type(ncvariable) :: var
90
91
92 ! construct filename
93 fname = trim(filelatlon(ii))
94
95 nc = ncdataset(fname, "r")
96
97 ! -------------------------------------------------------------------------
98 ! READ LEVEL LATITUDE / LONGITUDE
99 ! -------------------------------------------------------------------------
100 var = nc%getVariable(trim(lat_var_name))
101 call var%getData(dummy)
102 ! consistency check
103 if ((size(dummy, dim = 1) .NE. level%nrows) .or. &
104 (size(dummy, dim = 2) .NE. level%ncols)) then
105 call error_message(' ***ERROR: subroutine mo_read_latlon: size mismatch in latlon file for ', trim(level_name), &
106 ' in domain ', trim(adjustl(num2str(ii))), '!', raise=.false.)
107 call error_message(' Latitude expected to have following dimensions ... rows:', &
108 trim(adjustl(num2str(level%nrows))), ', cols:', trim(adjustl(num2str(level%ncols))), raise=.false.)
109 call error_message(' Latitude provided ... rows:', &
110 trim(adjustl(num2str(size(dummy, dim = 1)))), ', cols:', trim(adjustl(num2str(size(dummy, dim = 2)))))
111 end if
112 level%y = dummy
113
114 var = nc%getVariable(trim(lon_var_name))
115 call var%getData(dummy)
116 ! consistency check
117 if ((size(dummy, dim = 1) .NE. level%nrows) .or. &
118 (size(dummy, dim = 2) .NE. level%ncols)) then
119 call error_message(' ***ERROR: subroutine mo_read_latlon: size mismatch in latlon file for ', trim(level_name), &
120 ' in domain ', trim(adjustl(num2str(ii))), '!', raise=.false.)
121 call error_message(' Longitude expected to have following dimensions ... rows:', &
122 trim(adjustl(num2str(level%nrows))), ', cols:', trim(adjustl(num2str(level%ncols))), raise=.false.)
123 call error_message(' Longitude provided ... rows:', &
124 trim(adjustl(num2str(size(dummy, dim = 1)))), ', cols:', trim(adjustl(num2str(size(dummy, dim = 2)))))
125 end if
126 level%x = dummy
127
128 call nc%close()
129
130 end subroutine read_latlon
131
132END MODULE mo_read_latlon
Provides common types needed by mHM, mRM and/or mpr.
Provides structures needed by mHM, mRM and/or mpr.
character(256), dimension(:), allocatable, public filelatlon
reading latitude and longitude coordinates for each domain
subroutine, public read_latlon(ii, lon_var_name, lat_var_name, level_name, level)
reads latitude and longitude coordinates