VTK  9.5.20250806
vtkGeoTransform.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright 2008 Sandia Corporation
3// SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-USGov
15#ifndef vtkGeoTransform_h
16#define vtkGeoTransform_h
17
19#include "vtkGeovisCoreModule.h" // For export macro
20
21VTK_ABI_NAMESPACE_BEGIN
23
36class VTKGEOVISCORE_EXPORT vtkGeoTransform : public vtkAbstractTransform
37{
38public:
40 void PrintSelf(ostream& os, vtkIndent indent) override;
42
44
50 void SetSourceProjection(const char* proj);
51 vtkGetObjectMacro(SourceProjection, vtkGeoProjection);
53
55
65 vtkSetMacro(TransformZCoordinate, bool);
66 vtkGetMacro(TransformZCoordinate, bool);
68
70
76 void SetDestinationProjection(const char* proj);
77 vtkGetObjectMacro(DestinationProjection, vtkGeoProjection);
79
83 void TransformPoints(vtkPoints* src, vtkPoints* dst) override;
84
91 vtkDataArray* dstNms, vtkDataArray* srcVrs, vtkDataArray* dstVrs, int nOptionalVectors = 0,
92 vtkDataArray** srcVrsArr = nullptr, vtkDataArray** dstVrsArr = nullptr) override;
93
97 void Inverse() override;
98
100
104 void InternalTransformPoint(const float in[3], float out[3]) override;
105 void InternalTransformPoint(const double in[3], double out[3]) override;
107
109
116 static int ComputeUTMZone(double lon, double lat);
117 static int ComputeUTMZone(double* lonlat) { return ComputeUTMZone(lonlat[0], lonlat[1]); }
119
121
128 const float in[3], float out[3], float derivative[3][3]) override;
130 const double in[3], double out[3], double derivative[3][3]) override;
132
137
138protected:
141
143 void InternalTransformPoints(double* ptsInOut, vtkIdType numPts, int stride);
147
148private:
149 vtkGeoTransform(const vtkGeoTransform&) = delete;
150 void operator=(const vtkGeoTransform&) = delete;
151};
152
153VTK_ABI_NAMESPACE_END
154#endif // vtkGeoTransform_h
superclass for all geometric transformations
abstract superclass for arrays of numeric data
Represent a projection from a sphere to a plane.
A transformation between two geographic coordinate systems.
void TransformPointsNormalsVectors(vtkPoints *src, vtkPoints *dst, vtkDataArray *srcNms, vtkDataArray *dstNms, vtkDataArray *srcVrs, vtkDataArray *dstVrs, int nOptionalVectors=0, vtkDataArray **srcVrsArr=nullptr, vtkDataArray **dstVrsArr=nullptr) override
Another method to transform points, normals, and vectors all at once.
void TransformPoints(vtkPoints *src, vtkPoints *dst) override
Transform many points at once.
void Inverse() override
Invert the transformation.
void InternalTransformDerivative(const double in[3], double out[3], double derivative[3][3]) override
This will transform a point and, at the same time, calculate a 3x3 Jacobian matrix that provides the ...
void SetSourceProjection(vtkGeoProjection *source)
The source geographic projection, which can be set using an external vtkGeoProjection,...
void InternalTransformPoint(const float in[3], float out[3]) override
This will calculate the transformation without calling Update.
void InternalDeepCopy(vtkAbstractTransform *) override
Perform any subclass-specific DeepCopy.
static vtkGeoTransform * New()
vtkAbstractTransform * MakeTransform() override
Make another transform of the same type.
void InternalTransformDerivative(const float in[3], float out[3], float derivative[3][3]) override
This will transform a point and, at the same time, calculate a 3x3 Jacobian matrix that provides the ...
void SetSourceProjection(const char *proj)
The source geographic projection, which can be set using an external vtkGeoProjection,...
~vtkGeoTransform() override
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void InternalTransformPoint(const double in[3], double out[3]) override
This will calculate the transformation without calling Update.
static int ComputeUTMZone(double *lonlat)
Computes Universal Transverse Mercator (UTM) zone given the longitude and latitude of the point.
vtkGeoProjection * SourceProjection
void SetDestinationProjection(vtkGeoProjection *dest)
The target geographic projection, which can be set using an external vtkGeoProjection,...
static int ComputeUTMZone(double lon, double lat)
Computes Universal Transverse Mercator (UTM) zone given the longitude and latitude of the point.
void SetDestinationProjection(const char *proj)
The target geographic projection, which can be set using an external vtkGeoProjection,...
vtkGeoProjection * DestinationProjection
void InternalTransformPoints(double *ptsInOut, vtkIdType numPts, int stride)
a simple class to control print indentation
Definition vtkIndent.h:108
represent and manipulate 3D points
Definition vtkPoints.h:139
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition vtkType.h:332