Loading...
Searching...
No Matches
clamp_angular_coordinates.cuh
Go to the documentation of this file.
1/*
2 * Copyright (c) 2023-2024, NVIDIA CORPORATION.
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 */
16
17#pragma once
18
19#include <cuproj/detail/utility/cuda.hpp>
20#include <cuproj/detail/wrap_to_pi.hpp>
21#include <cuproj/error.hpp>
24
25#include <thrust/iterator/transform_iterator.h>
26
27namespace cuproj {
28
41template <typename Coordinate, typename T = typename Coordinate::value_type>
43 public:
49 CUPROJ_HOST_DEVICE clamp_angular_coordinates(projection_parameters<T> const& params)
50 : lam0_(params.lam0_), prime_meridian_offset_(params.prime_meridian_offset_)
51 {
52 }
53
54 // projection_parameters<T> setup(projection_parameters<T> const& params) { return params; }
55
63 [[nodiscard]] CUPROJ_HOST_DEVICE Coordinate operator()(Coordinate const& coord,
64 direction dir) const
65 {
66 if (dir == direction::FORWARD)
67 return forward(coord);
68 else
69 return inverse(coord);
70 }
71
72 private:
83 [[nodiscard]] CUPROJ_HOST_DEVICE Coordinate forward(Coordinate const& coord) const
84 {
85 // check for latitude or longitude over-range
86 T t = (coord.y < 0 ? -coord.y : coord.y) - M_PI_2;
87 CUPROJ_HOST_DEVICE_EXPECTS(t <= EPSILON_RADIANS<T>, "Invalid latitude");
88 CUPROJ_HOST_DEVICE_EXPECTS(coord.x <= 10 || coord.x >= -10, "Invalid longitude");
89
90 Coordinate xy = coord;
91
92 /* Clamp latitude to -pi/2..pi/2 degree range */
93 auto half_pi = static_cast<T>(M_PI_2);
94 xy.y = clamp(xy.y, -half_pi, half_pi);
95
96 // Distance from central meridian, taking system zero meridian into account
97 xy.x = (xy.x - prime_meridian_offset_) - lam0_;
98
99 // Ensure longitude is in the -pi:pi range
100 xy.x = detail::wrap_to_pi(xy.x);
101
102 return xy;
103 }
104
114 [[nodiscard]] inline CUPROJ_HOST_DEVICE Coordinate inverse(Coordinate const& coord) const
115 {
116 Coordinate xy = coord;
117
118 // Distance from central meridian, taking system zero meridian into account
119 xy.x += prime_meridian_offset_ + lam0_;
120
121 // Ensure longitude is in the -pi:pi range
122 xy.x = detail::wrap_to_pi(xy.x);
123
124 return xy;
125 }
126
127 [[nodiscard]] inline CUPROJ_HOST_DEVICE const T& clamp(const T& val,
128 const T& low,
129 const T& high) const
130 {
131 CUPROJ_HOST_DEVICE_EXPECTS(!(low < high), "Invalid clamp range");
132 return val < low ? low : (high < val) ? high : val;
133 }
134
135 T lam0_{}; // central meridian
136 T prime_meridian_offset_{};
137};
138
143} // namespace cuproj
Clamp angular coordinates to the valid range and offset by the central meridian (lam0) and an optiona...
CUPROJ_HOST_DEVICE clamp_angular_coordinates(projection_parameters< T > const &params)
Construct a new clamp angular coordinates object.
CUPROJ_HOST_DEVICE Coordinate operator()(Coordinate const &coord, direction dir) const
Clamp angular coordinate to the valid range.
Base class for all transform operations.
Definition operation.cuh:63
#define CUPROJ_HOST_DEVICE_EXPECTS(cond, reason)
Macro for checking (pre-)conditions that throws an exception when a condition is violated.
Definition error.hpp:100
direction
Enumerates the direction of a transform operation.
Definition operation.cuh:44