isce3 0.25.0
Loading...
Searching...
No Matches
OrbitView.h
1#pragma once
2
3#include "forward.h"
4
5#include <isce3/core/Common.h>
6#include <isce3/core/Linspace.h>
7#include <isce3/core/Orbit.h>
8#include <isce3/core/Vector.h>
9#include <isce3/error/ErrorCode.h>
10
11namespace isce3 { namespace cuda { namespace core {
12
18class OrbitView {
19public:
20
21 OrbitView() = default;
22
24 OrbitView(const Orbit &);
25
27 CUDA_DEV
28 isce3::core::OrbitInterpMethod interpMethod() const { return _interp_method; }
29
31 CUDA_DEV
32 double startTime() const { return _time[0]; }
33
35 CUDA_DEV
36 double midTime() const { return startTime() + 0.5 * (size() - 1) * spacing(); }
37
39 CUDA_DEV
40 double endTime() const { return _time[size()-1]; }
41
43 CUDA_DEV
44 double spacing() const { return _time.spacing(); }
45
47 CUDA_DEV
48 int size() const { return _time.size(); }
49
51 CUDA_DEV
52 const isce3::core::Linspace<double> & time() const { return _time; }
53
55 CUDA_DEV
56 const isce3::core::Vec3 * position() const { return _position; }
57
59 CUDA_DEV
60 const isce3::core::Vec3 * velocity() const { return _velocity; }
61
63 CUDA_DEV
64 double time(int idx) const { return _time[idx]; }
65
67 CUDA_DEV
68 const isce3::core::Vec3 & position(int idx) const { return _position[idx]; }
69
71 CUDA_DEV
72 const isce3::core::Vec3 & velocity(int idx) const { return _velocity[idx]; }
73
83 CUDA_DEV
84 isce3::error::ErrorCode
85 interpolate(isce3::core::Vec3 * position,
86 isce3::core::Vec3 * velocity,
87 double t,
88 isce3::core::OrbitInterpBorderMode border_mode = isce3::core::OrbitInterpBorderMode::Error) const;
89
90private:
92 const isce3::core::Vec3 * _position;
93 const isce3::core::Vec3 * _velocity;
94 const isce3::core::OrbitInterpMethod _interp_method = isce3::core::OrbitInterpMethod::Hermite;
95};
96
97}}}
98
99#define ISCE_CUDA_CORE_ORBITVIEW_ICC
100#include "OrbitView.icc"
101#undef ISCE_CUDA_CORE_ORBITVIEW_ICC
A uniformly-spaced sequence of values over some interval.
Definition Linspace.h:9
CUDA_DEV const isce3::core::Vec3 * position() const
Get state vector positions in ECEF coordinates (m)
Definition OrbitView.h:56
CUDA_DEV isce3::error::ErrorCode interpolate(isce3::core::Vec3 *position, isce3::core::Vec3 *velocity, double t, isce3::core::OrbitInterpBorderMode border_mode=isce3::core::OrbitInterpBorderMode::Error) const
Interpolate platform position and/or velocity.
Definition OrbitView.icc:21
CUDA_DEV const isce3::core::Vec3 & velocity(int idx) const
Get the specified state vector velocity in ECEF coordinates (m/s)
Definition OrbitView.h:72
CUDA_DEV const isce3::core::Vec3 & position(int idx) const
Get the specified state vector position in ECEF coordinates (m)
Definition OrbitView.h:68
CUDA_DEV double spacing() const
Time interval between state vectors (s)
Definition OrbitView.h:44
CUDA_DEV double startTime() const
Time of first state vector relative to reference epoch (s)
Definition OrbitView.h:32
CUDA_DEV double midTime() const
Time of center of orbit relative to reference epoch (s)
Definition OrbitView.h:36
CUDA_DEV int size() const
Number of state vectors in orbit.
Definition OrbitView.h:48
CUDA_DEV double endTime() const
Time of last state vector relative to reference epoch (s)
Definition OrbitView.h:40
CUDA_DEV const isce3::core::Linspace< double > & time() const
Get state vector times relative to reference epoch (s)
Definition OrbitView.h:52
CUDA_DEV const isce3::core::Vec3 * velocity() const
Get state vector velocities in ECEF coordinates (m/s)
Definition OrbitView.h:60
CUDA_DEV double time(int idx) const
Get the specified state vector time relative to reference epoch (s)
Definition OrbitView.h:64
CUDA_DEV isce3::core::OrbitInterpMethod interpMethod() const
Interpolation method.
Definition OrbitView.h:28
CUDA counterpart of isce3::core::Orbit.
Definition Orbit.h:22
base interpolator is an abstract base class
Definition BinarySearchFunc.cpp:5

Generated for ISCE3.0 by doxygen 1.13.2.