libcuproj  24.02.00
Loading...
Searching...
No Matches
projection.cuh
Go to the documentation of this file.
1/*
2 * Copyright (c) 2023, 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/pipeline.cuh>
20#include <cuproj/ellipsoid.hpp>
23
24#include <rmm/cuda_stream_view.hpp>
25#include <rmm/exec_policy.hpp>
26
27#include <thrust/device_vector.h>
28#include <thrust/transform.h>
29
30#include <iterator>
31#include <type_traits>
32
33namespace cuproj {
34
50template <typename Coordinate, typename T = typename Coordinate::value_type>
52 public:
60 projection(std::vector<operation_type> const& operations,
61 projection_parameters<T> const& params,
62 direction dir = direction::FORWARD)
63 : params_(params), constructed_direction_(dir)
64 {
65 setup(operations);
66 }
67
79 template <class InputCoordIter, class OutputCoordIter>
80 void transform(InputCoordIter first,
81 InputCoordIter last,
82 OutputCoordIter result,
83 direction dir,
84 rmm::cuda_stream_view stream = rmm::cuda_stream_default) const
85 {
86 dir = (constructed_direction_ == direction::FORWARD) ? dir : reverse(dir);
87
88 if (dir == direction::FORWARD) {
89 auto pipe = detail::pipeline<Coordinate, direction::FORWARD>{
90 params_, operations_.data().get(), operations_.size()};
91 thrust::transform(rmm::exec_policy(stream), first, last, result, pipe);
92 } else {
93 auto pipe = detail::pipeline<Coordinate, direction::INVERSE>{
94 params_, operations_.data().get(), operations_.size()};
95 thrust::transform(rmm::exec_policy(stream), first, last, result, pipe);
96 }
97 }
98
99 private:
100 void setup(std::vector<operation_type> const& operations)
101 {
102 std::for_each(operations.begin(), operations.end(), [&](auto const& op) {
103 switch (op) {
104 case operation_type::TRANSVERSE_MERCATOR: {
105 auto op = transverse_mercator<Coordinate>{params_};
106 params_ = op.setup(params_);
107 break;
108 }
109 // TODO: some ops don't have setup. Should we make them all have setup?
110 default: break;
111 }
112 });
113
114 operations_.resize(operations.size());
115 thrust::copy(operations.begin(), operations.end(), operations_.begin());
116 }
117
118 thrust::device_vector<operation_type> operations_;
119 projection_parameters<T> params_;
120 direction constructed_direction_{direction::FORWARD};
121};
122
127} // namespace cuproj
A projection transforms coordinates between coordinate reference systems.
projection(std::vector< operation_type > const &operations, projection_parameters< T > const &params, direction dir=direction::FORWARD)
Construct a new projection object.
void transform(InputCoordIter first, InputCoordIter last, OutputCoordIter result, direction dir, rmm::cuda_stream_view stream=rmm::cuda_stream_default) const
Transform a range of coordinates.
direction
Enumerates the direction of a transform operation.
Definition operation.cuh:44
direction reverse(direction dir)
Returns the opposite of a direction.
Definition operation.cuh:47