Isis 3 Programmer Reference
CSMCamera.cpp
1
7/* SPDX-License-Identifier: CC0-1.0 */
8
9#include "CSMCamera.h"
10#include "CSMSkyMap.h"
11
12#include <fstream>
13#include <iostream>
14#include <iomanip>
15
16#include <QDebug>
17#include <QList>
18#include <QPointF>
19#include <QString>
20
21#include "Affine.h"
22#include "Blob.h"
23#include "Constants.h"
24#include "Displacement.h"
25#include "Distance.h"
26#include "FileName.h"
27#include "IException.h"
28#include "IString.h"
29#include "iTime.h"
30#include "Latitude.h"
31#include "Longitude.h"
32#include "LinearAlgebra.h"
33#include "NaifStatus.h"
34#include "SpecialPixel.h"
35#include "SurfacePoint.h"
36#include "Table.h"
37#include "SensorUtilities.h"
38
39#include "csm/Warning.h"
40#include "csm/Error.h"
41#include "csm/Plugin.h"
42#include "csm/Ellipsoid.h"
43#include "csm/SettableEllipsoid.h"
44
45#include <nlohmann/json.hpp>
46using json = nlohmann::json;
47
48using namespace std;
49
50namespace Isis {
51
52json stateAsJson(std::string modelState);
53void sanitize(std::string &input);
54
63 Blob state("CSMState", "String");
64 cube.read(state);
65 PvlObject &blobLabel = state.Label();
66
67 QString pluginName = blobLabel.findKeyword("PluginName")[0];
68 QString modelName = blobLabel.findKeyword("ModelName")[0];
69 QString stateString = QString::fromUtf8(state.getBuffer(), state.Size());
70 init(cube, pluginName, modelName, stateString);
71 }
72
73
82 void CSMCamera::init(Cube &cube, QString pluginName, QString modelName, QString stateString) {
83 const csm::Plugin *plugin = csm::Plugin::findPlugin(pluginName.toStdString());
84 if (!plugin) {
85 QStringList availablePlugins;
86 for (const csm::Plugin *plugin: csm::Plugin::getList()) {
87 availablePlugins.append(QString::fromStdString(plugin->getPluginName()));
88 }
89 QString msg = "Failed to find plugin [" + pluginName + "] for image [" + cube.fileName() +
90 "]. Check that the corresponding CSM plugin library is in the directory "
91 "specified by your IsisPreferences. Loaded plugins [" +
92 availablePlugins.join(", ") + "].";
93 throw IException(IException::User, msg, _FILEINFO_);
94 }
95 if (!plugin->canModelBeConstructedFromState(modelName.toStdString(), stateString.toStdString())) {
96 QString msg = "CSM state string attached to image [" + cube.fileName() + "] cannot "
97 "be converted to a [" + modelName + "] using [" + pluginName + "].";
98 throw IException(IException::Programmer, msg, _FILEINFO_);
99 }
100 m_model = dynamic_cast<csm::RasterGM*>(plugin->constructModelFromState(stateString.toStdString()));
101 // If the dynamic cast failed, raise an exception
102 if (!m_model) {
103 QString msg = "Failed to convert CSM Model to RasterGM.";
104 throw IException(IException::Programmer, msg, _FILEINFO_);
105 }
106
107 m_instrumentNameLong = QString::fromStdString(m_model->getSensorIdentifier());
108 m_instrumentNameShort = QString::fromStdString(m_model->getSensorIdentifier());
109 m_spacecraftNameLong = QString::fromStdString(m_model->getPlatformIdentifier());
110 m_spacecraftNameShort = QString::fromStdString(m_model->getPlatformIdentifier());
111
112 QString timeString = QString::fromStdString(m_model->getReferenceDateAndTime());
113 // Strip the UTC time zone indicator for iTime
114 timeString.remove("Z");
115 m_refTime.setUtc(timeString);
116
117 setTarget(*cube.label());
118 }
119
120
130 bool CSMCamera::SetImage(const double sample, const double line) {
131 // Save off the line & sample
132 p_childSample = sample;
133 p_childLine = line;
134
135 csm::ImageCoord imagePt;
137 double achievedPrecision = 0;
138 csm::WarningList warnings;
139 csm::EcefLocus imageLocus;
140 try {
141 imageLocus = m_model->imageToRemoteImagingLocus(imagePt,
142 0.001,
143 &achievedPrecision,
144 &warnings);
145 }
146 catch (csm::Error &e) {
147 return false;
148 }
149
150 // Check for issues on the CSM end
151 if (achievedPrecision > 0.001) {
152 return false;
153 }
154 if (!warnings.empty()) {
155 for (csm::Warning warning : warnings) {
156 if (warning.getWarning() == csm::Warning::IMAGE_COORD_OUT_OF_BOUNDS){
157 return false;
158 }
159 }
160 }
161
162 // ISIS sensors work in Kilometers, CSM works in meters
163 std::vector<double> obsPosition = {imageLocus.point.x / 1000.0,
164 imageLocus.point.y / 1000.0,
165 imageLocus.point.z / 1000.0};
166 std::vector<double> locusVec = {imageLocus.direction.x,
167 imageLocus.direction.y,
168 imageLocus.direction.z};
169
170 // Save off the look vector
171 SetLookDirection(locusVec);
172 if (!m_et) {
173 m_et = new iTime();
174 }
175 *m_et = m_refTime + m_model->getImageTime(imagePt);
176 if (target()->isSky()) {
177 target()->shape()->setHasIntersection(false);
178 return true;
179 }
180
181 // Check for a ground intersection
182 if(!target()->shape()->intersectSurface(obsPosition, locusVec)) {
183 return false;
184 }
185
186 p_pointComputed = true;
187 return true;
188 }
189
190
201 bool CSMCamera::SetUniversalGround(const double latitude, const double longitude) {
202 return SetGround(
203 Latitude(latitude, Angle::Degrees),
204 Longitude(longitude, Angle::Degrees));
205 }
206
207
218 bool CSMCamera::SetUniversalGround(const double latitude, const double longitude, double radius) {
219 return SetGround(SurfacePoint(
220 Latitude(latitude, Angle::Degrees),
221 Longitude(longitude, Angle::Degrees),
222 Distance(radius, Distance::Meters)));
223 }
224
225
236 bool CSMCamera::SetGround(Latitude latitude, Longitude longitude) {
237 ShapeModel *shape = target()->shape();
238 Distance localRadius;
239
240 if (shape->name() != "Plane") { // this is the normal behavior
241 localRadius = LocalRadius(latitude, longitude);
242 }
243 else {
244 localRadius = Distance(latitude.degrees(),Distance::Kilometers);
245 latitude = Latitude(0.,Angle::Degrees);
246 }
247
248 if (!localRadius.isValid()) {
249 target()->shape()->clearSurfacePoint();
250 return false;
251 }
252
253 return SetGround(SurfacePoint(latitude, longitude, localRadius));
254 }
255
256
265 bool CSMCamera::SetRightAscensionDeclination(const double ra, const double dec) {
266 double raRad = ra * DEG2RAD;
267 double decRad = dec * DEG2RAD;
268
269 // Make the radius bigger, some multiple of the body radius -or- use sensor position at the reference point
270 SensorUtilities::GroundPt3D sphericalPt = {decRad, raRad, 1};
271 SensorUtilities::Vec rectPt = SensorUtilities::sphericalToRect(sphericalPt);
272
273 csm::EcefCoord lookPt = {rectPt.x, rectPt.y, rectPt.z};
274 double achievedPrecision = 0;
275 csm::WarningList warnings;
276 bool validBackProject;
277 csm::ImageCoord imagePt;
278
279 try {
280 imagePt = m_model->groundToImage(lookPt, 0.01, &achievedPrecision, &warnings);
281 double sample;
282 double line;
283 csmToIsisPixel(imagePt, line, sample);
284 validBackProject = SetImage(sample, line);
285 }
286 catch (csm::Error &e) {
287 validBackProject = false;
288 }
289
290 return validBackProject;
291 }
292
293
306 bool CSMCamera::SetLookDirection(const std::vector<double> lookB) {
307 // The look vector must be in the camera coordinate system
308
309 // This memcpy does:
310 // m_lookB[0] = lookB[0];
311 // m_lookB[1] = lookB[1];
312 // m_lookB[2] = lookB[2];
313 memcpy(m_lookB, &lookB[0], sizeof(double) * 3);
314 m_newLookB = true;
315
316 // Don't try to intersect the sky
317 if (target()->isSky()) {
318 target()->shape()->setHasIntersection(false);
319 return false;
320 }
321
322 return false;
323 }
324
325
337 double orgLine = Line();
338 double orgSample = Sample();
339 double orgDec = Declination();
340 double orgRa = RightAscension();
341
342 SetRightAscensionDeclination(orgRa, orgDec + (2 * RaDecResolution()));
343 double y = Line() - orgLine;
344 double x = Sample() - orgSample;
345 double celestialNorthClockAngle = atan2(-y, x) * 180.0 / Isis::PI;
346 celestialNorthClockAngle = 90.0 - celestialNorthClockAngle;
347
348 if (celestialNorthClockAngle < 0.0) {
349 celestialNorthClockAngle += 360.0;
350 }
351
352 SetImage(orgSample, orgLine);
353 return celestialNorthClockAngle;
354 }
355
356
365 bool CSMCamera::SetGround(const SurfacePoint & surfacePt) {
366 ShapeModel *shape = target()->shape();
367 if (!surfacePt.Valid()) {
368 shape->clearSurfacePoint();
369 return false;
370 }
371
372 bool validBackProject = true;
373
374 // Back project through the CSM model
375 csm::ImageCoord imagePt;
376 double achievedPrecision = 0;
377 csm::WarningList warnings;
378 csm::EcefCoord groundPt = isisToCsmGround(surfacePt);
379 try {
380 imagePt = m_model->groundToImage(groundPt, 0.01, &achievedPrecision, &warnings);
381 }
382 catch (csm::Error &e) {
383 validBackProject = false;
384 }
385 if (achievedPrecision > 0.01) {
386 validBackProject = false;
387 }
388 if (!warnings.empty()) {
389 for (csm::Warning warning : warnings) {
390 if (warning.getWarning() == csm::Warning::IMAGE_COORD_OUT_OF_BOUNDS){
391 validBackProject = false;
392 }
393 }
394 }
395
396 // Check for occlusion
397 double line, sample;
398 csmToIsisPixel(imagePt, line, sample);
399 csm::EcefLocus imageLocus = m_model->imageToRemoteImagingLocus(imagePt);
400 std::vector<double> sensorPosition = {imageLocus.point.x, imageLocus.point.y, imageLocus.point.z};
401 shape->clearSurfacePoint();
402 shape->intersectSurface(surfacePt,
403 sensorPosition,
404 true);
405 if (!shape->hasIntersection()) {
406 validBackProject = false;
407 }
408
409 // If the back projection was successful, then save it
410 if (validBackProject) {
411 m_lookB[0] = imageLocus.direction.x;
412 m_lookB[1] = imageLocus.direction.y;
413 m_lookB[2] = imageLocus.direction.z;
414 m_newLookB = true;
417 p_pointComputed = true;
418 shape->setHasIntersection(true);
419 if (!m_et) {
420 m_et = new iTime();
421 }
422 *m_et = m_refTime + m_model->getImageTime(imagePt);
423 return true;
424 }
425
426 // Otherwise reset
427 shape->clearSurfacePoint();
428 return false;
429 }
430
431
447 vector<double> imagePartials = ImagePartials();
448 return sqrt(imagePartials[0]*imagePartials[0] +
449 imagePartials[2]*imagePartials[2] +
450 imagePartials[4]*imagePartials[4]);
451 }
452
453 return Isis::Null;
454 }
455
456
468 vector<double> imagePartials = ImagePartials();
469 return sqrt(imagePartials[1] * imagePartials[1] +
470 imagePartials[3] * imagePartials[3] +
471 imagePartials[5] * imagePartials[5]);
472 }
473
474 return Isis::Null;
475 }
476
477
489 // Redo the line and sample resolution calculations because it avoids
490 // a call to ImagePartials which could be a costly call
491 vector<double> imagePartials = ImagePartials();
492 double lineRes = sqrt(imagePartials[0]*imagePartials[0] +
493 imagePartials[2]*imagePartials[2] +
494 imagePartials[4]*imagePartials[4]);
495 double sampRes = sqrt(imagePartials[1]*imagePartials[1] +
496 imagePartials[3]*imagePartials[3] +
497 imagePartials[5]*imagePartials[5]);
498 return (sampRes + lineRes) / 2.0;
499 }
500
501 return Isis::Null;
502 }
503
504
514 double CSMCamera::ObliqueLineResolution(bool useLocal) {
515 // CSM resolution is always the oblique resolution so just return it
516 return LineResolution();
517 }
518
519
530 // CSM resolution is always the oblique resolution so just return it
531 return SampleResolution();
532 }
533
534
545 // CSM resolution is always the oblique resolution so just return it
546 return DetectorResolution();
547 }
548
549
557 double CSMCamera::parentLine() const {
558 return p_alphaCube->AlphaLine(Line());
559 }
560
561
569 double CSMCamera::parentSample() const {
570 return p_alphaCube->AlphaSample(Sample());
571 }
572
573
582 std::vector<double> position = sensorPositionBodyFixed();
583 p[0] = position[0];
584 p[1] = position[1];
585 p[2] = position[2];
586 }
587
588
595 std::vector<double> CSMCamera::sensorPositionBodyFixed() const {
597 }
598
599
609 std::vector<double> CSMCamera::sensorPositionBodyFixed(double line, double sample) const {
610 csm::ImageCoord imagePt;
611 isisToCsmPixel(line, sample, imagePt);
612 csm::EcefCoord sensorPosition = m_model->getSensorPosition(imagePt);
613 // CSM uses meters, but ISIS wants this in Km
614 std::vector<double> result {
615 sensorPosition.x / 1000.0,
616 sensorPosition.y / 1000.0,
617 sensorPosition.z / 1000.0};
618 return result;
619 }
620
621
630 void CSMCamera::subSpacecraftPoint(double &lat, double &lon) {
632 }
633
634
645 void CSMCamera::subSpacecraftPoint(double &lat, double &lon, double line, double sample) {
646 // Get s/c position from CSM because it is vector from center of body to that
647 vector<double> sensorPosition = sensorPositionBodyFixed(line, sample);
648 SurfacePoint surfacePoint(
649 Displacement(sensorPosition[0], Displacement::Kilometers),
650 Displacement(sensorPosition[1], Displacement::Kilometers),
651 Displacement(sensorPosition[2], Displacement::Kilometers));
652 lat = surfacePoint.GetLatitude().degrees();
653 lon = surfacePoint.GetLongitude().degrees();
654 }
655
656
672 vector<double> CSMCamera::ImagePartials() {
674 }
675
676
699 vector<double> CSMCamera::ImagePartials(SurfacePoint groundPoint) {
700 csm::EcefCoord groundCoord = isisToCsmGround(groundPoint);
701 vector<double> groundPartials = m_model->computeGroundPartials(groundCoord);
702
703 // Jacobian format is
704 // line WRT X line WRT Y line WRT Z
705 // samp WRT X samp WRT Y samp WRT Z
706 LinearAlgebra::Matrix groundMatrix(2, 3);
707 groundMatrix(0,0) = groundPartials[0];
708 groundMatrix(0,1) = groundPartials[1];
709 groundMatrix(0,2) = groundPartials[2];
710 groundMatrix(1,0) = groundPartials[3];
711 groundMatrix(1,1) = groundPartials[4];
712 groundMatrix(1,2) = groundPartials[5];
713
714 LinearAlgebra::Matrix imageMatrix = LinearAlgebra::pseudoinverse(groundMatrix);
715
716 vector<double> imagePartials = {imageMatrix(0,0),
717 imageMatrix(0,1),
718 imageMatrix(1,0),
719 imageMatrix(1,1),
720 imageMatrix(2,0),
721 imageMatrix(2,1)};
722 return imagePartials;
723 }
724
725
742 vector<double> CSMCamera::GroundPartials() {
744 }
745
746
765 vector<double> CSMCamera::GroundPartials(SurfacePoint groundPoint) {
766 csm::EcefCoord groundCoord = isisToCsmGround(groundPoint);
767 vector<double> groundPartials = m_model->computeGroundPartials(groundCoord);
768 return groundPartials;
769 }
770
771
778 Target *target = new Target(label);
779
780 // get radii from CSM
781 csm::Ellipsoid targetEllipsoid = csm::SettableEllipsoid::getEllipsoid(m_model);
782 std::vector<Distance> radii = {Distance(targetEllipsoid.getSemiMajorRadius(), Distance::Meters),
783 Distance(targetEllipsoid.getSemiMajorRadius(), Distance::Meters),
784 Distance(targetEllipsoid.getSemiMinorRadius(), Distance::Meters)};
785 target->setRadii(radii);
786
787 // Target needs to be able to access the camera to do things like
788 // compute resolution
789 target->setSpice(this);
790
791 if (m_target) {
792 delete m_target;
793 m_target = nullptr;
794 }
795
797 }
798
799
809 void CSMCamera::isisToCsmPixel(double line, double sample, csm::ImageCoord &csmPixel) const {
810 csmPixel.line = line - 0.5;
811 csmPixel.samp = sample - 0.5;
812 }
813
814
824 void CSMCamera::csmToIsisPixel(csm::ImageCoord csmPixel, double &line, double &sample) const {
825 line = csmPixel.line + 0.5;
826 sample = csmPixel.samp + 0.5;
827 }
828
829
840 csm::EcefCoord CSMCamera::isisToCsmGround(const SurfacePoint &groundPt) const {
841 return csm::EcefCoord(groundPt.GetX().meters(),
842 groundPt.GetY().meters(),
843 groundPt.GetZ().meters());
844 }
845
846
857 SurfacePoint CSMCamera::csmToIsisGround(const csm::EcefCoord &groundPt) const {
861 }
862
863
869 double CSMCamera::PhaseAngle() const {
870 csm::EcefCoord groundPt = isisToCsmGround(GetSurfacePoint());
871 csm::EcefVector sunEcefVec = m_model->getIlluminationDirection(groundPt);
872 // ISIS wants the position of the sun, not just the vector from the ground
873 // point to the sun. So, we approximate this by adding in the ground point.
874 // ISIS wants this in Km so convert
875 std::vector<double> sunVec = {
876 (groundPt.x - sunEcefVec.x) / 1000.0,
877 (groundPt.y - sunEcefVec.y) / 1000.0,
878 (groundPt.z - sunEcefVec.z) / 1000.0};
879 return target()->shape()->phaseAngle(sensorPositionBodyFixed(), sunVec);
880 }
881
882
889 return target()->shape()->emissionAngle(sensorPositionBodyFixed());
890 }
891
892
899 csm::EcefCoord groundPt = isisToCsmGround(GetSurfacePoint());
900 csm::EcefVector sunEcefVec = m_model->getIlluminationDirection(groundPt);
901 // ISIS wants the position of the sun, not just the vector from the ground
902 // point to the sun. So, we approximate this by adding in the ground point.
903 // ISIS wants this in Km so convert
904 std::vector<double> sunVec = {
905 (groundPt.x - sunEcefVec.x) / 1000.0,
906 (groundPt.y - sunEcefVec.y) / 1000.0,
907 (groundPt.z - sunEcefVec.z) / 1000.0};
908 return target()->shape()->incidenceAngle(sunVec);
909 }
910
911
919 std::vector<double> sensorPosition = sensorPositionBodyFixed();
920 SurfacePoint groundPoint = GetSurfacePoint();
921
922 std::vector<double> sensorToGround = {
923 groundPoint.GetX().kilometers() - (sensorPosition[0]),
924 groundPoint.GetY().kilometers() - (sensorPosition[1]),
925 groundPoint.GetZ().kilometers() - (sensorPosition[2])};
926
927 return sqrt(
928 sensorToGround[0] * sensorToGround[0] +
929 sensorToGround[1] * sensorToGround[1] +
930 sensorToGround[2] * sensorToGround[2]);
931 }
932
933
941 std::vector<double> sensorPosition = sensorPositionBodyFixed();
942 return sqrt(
943 sensorPosition[0] * sensorPosition[0] +
944 sensorPosition[1] * sensorPosition[1] +
945 sensorPosition[2] * sensorPosition[2]);
946 }
947
948
956 std::vector<int> CSMCamera::getParameterIndices(csm::param::Set paramSet) const {
957 return m_model->getParameterSetIndices(paramSet);
958 }
959
960
968 std::vector<int> CSMCamera::getParameterIndices(csm::param::Type paramType) const {
969 std::vector<int> parameterIndices;
970 for (int i = 0; i < m_model->getNumParameters(); i++) {
971 if (m_model->getParameterType(i) == paramType) {
972 parameterIndices.push_back(i);
973 }
974 }
975 return parameterIndices;
976 }
977
978
986 std::vector<int> CSMCamera::getParameterIndices(QStringList paramList) const {
987 std::vector<int> parameterIndices;
988 QStringList failedParams;
989 for (int i = 0; i < paramList.size(); i++) {
990 bool found = false;
991 for (int j = 0; j < m_model->getNumParameters(); j++) {
992 if (QString::compare(QString::fromStdString(m_model->getParameterName(j)).trimmed(),
993 paramList[i].trimmed(),
994 Qt::CaseInsensitive) == 0) {
995 parameterIndices.push_back(j);
996 found = true;
997 break;
998 }
999 }
1000
1001 if (!found) {
1002 failedParams.push_back(paramList[i]);
1003 }
1004 }
1005
1006 if (!failedParams.empty()) {
1007 QString msg = "Failed to find indices for the following parameters [" +
1008 failedParams.join(",") + "].";
1009 throw IException(IException::User, msg, _FILEINFO_);
1010 }
1011 return parameterIndices;
1012 }
1013
1014
1021 void CSMCamera::applyParameterCorrection(int index, double correction) {
1022 double currentValue = m_model->getParameterValue(index);
1023 m_model->setParameterValue(index, currentValue + correction);
1024 }
1025
1026
1033 double CSMCamera::getParameterCovariance(int index1, int index2) {
1034 return m_model->getParameterCovariance(index1, index2);
1035 }
1036
1037
1038 vector<double> CSMCamera::getSensorPartials(int index, SurfacePoint groundPoint) {
1039 // csm::SensorPartials holds (line, sample) in order for each parameter
1040 csm::EcefCoord groundCoord = isisToCsmGround(groundPoint);
1041 std::pair<double, double> partials = m_model->computeSensorPartials(index, groundCoord);
1042 vector<double> partialsVector = {partials.first, partials.second};
1043
1044 return partialsVector;
1045 }
1046
1047
1055 QString CSMCamera::getParameterName(int index) {
1056 return QString::fromStdString(m_model->getParameterName(index));
1057 }
1058
1059
1068 return m_model->getParameterValue(index);
1069 }
1070
1071
1079 QString CSMCamera::getParameterUnits(int index) {
1080 return QString::fromStdString(m_model->getParameterUnits(index));
1081 }
1082
1083
1090 return QString::fromStdString(m_model->getModelState());
1091 }
1092
1093
1102 void CSMCamera::setTime(const iTime &time) {
1103 QString msg = "Setting the image time is not supported for CSM camera models";
1104 throw IException(IException::Programmer, msg, _FILEINFO_);
1105 }
1106
1107
1118 void CSMCamera::subSolarPoint(double &lat, double &lon) {
1119 QString msg = "Sub solar point is not supported for CSM camera models";
1120 throw IException(IException::Programmer, msg, _FILEINFO_);
1121 }
1122
1123
1132 QString msg = "Pixel Field of View is not supported for CSM camera models";
1133 throw IException(IException::User, msg, _FILEINFO_);
1134 }
1135
1136
1145 void CSMCamera::sunPosition(double p[3]) const {
1146 QString msg = "Sun position is not supported for CSM camera models";
1147 throw IException(IException::Programmer, msg, _FILEINFO_);
1148 }
1149
1150
1160 QString msg = "Sun position is not supported for CSM camera models";
1161 throw IException(IException::Programmer, msg, _FILEINFO_);
1162 }
1163
1164
1175 QString msg = "Instrument position is not supported for CSM camera models";
1176 throw IException(IException::Programmer, msg, _FILEINFO_);
1177 }
1178
1189 QString msg = "Body orientation is not supported for CSM camera models";
1190 throw IException(IException::Programmer, msg, _FILEINFO_);
1191 }
1192
1193
1204 QString msg = "Instrument orientation is not supported for CSM camera models";
1205 throw IException(IException::Programmer, msg, _FILEINFO_);
1206 }
1207
1208
1218 QString msg = "Solar longitude is not supported for CSM camera models";
1219 throw IException(IException::Programmer, msg, _FILEINFO_);
1220 }
1221
1222
1231 QString msg = "Solar distance is not supported for CSM camera models";
1232 throw IException(IException::Programmer, msg, _FILEINFO_);
1233 }
1234
1235
1242 double precision;
1243 csm::EcefLocus locus = m_model->imageToRemoteImagingLocus(csm::ImageCoord(p_childLine, p_childSample), 0.00001, &precision);
1244 csm::EcefVector v = locus.direction;
1245 SensorUtilities::GroundPt3D sphere_v = SensorUtilities::rectToSpherical({v.x, v.y, v.z});
1246 double lon = sphere_v.lon;
1247 if (lon < 0) {
1248 lon += 2 * PI;
1249 }
1250 return lon * RAD2DEG;
1251 }
1252
1253
1260 double precision;
1261 csm::EcefLocus locus = m_model->imageToRemoteImagingLocus(csm::ImageCoord(p_childLine, p_childSample), 0.00001, &precision);
1262 csm::EcefVector v = locus.direction;
1263 SensorUtilities::GroundPt3D sphere_v = SensorUtilities::rectToSpherical({v.x, v.y, v.z});
1264 return sphere_v.lat * RAD2DEG;
1265 }
1266
1267 json stateAsJson(std::string modelState) {
1268 // Remove special characters from string
1269 sanitize(modelState);
1270
1271 std::size_t foundFirst = modelState.find_first_of("{");
1272 std::size_t foundLast = modelState.find_last_of("}");
1273
1274 if (foundFirst == std::string::npos) {
1275 foundFirst = 0;
1276 }
1277 return json::parse(modelState.begin() + foundFirst, modelState.begin() + foundLast + 1);
1278 }
1279
1280 void sanitize(std::string &input){
1281 // Replaces characters from the string that are not printable with newlines
1282 std::replace_if(input.begin(), input.end(), [](int c){return !::isprint(c);}, '\n');
1283 }
1284}
double AlphaSample(double betaSample)
Returns an alpha sample given a beta sample.
Definition AlphaCube.h:121
double BetaLine(double alphaLine)
Returns a beta line given an alpha line.
Definition AlphaCube.h:133
double AlphaLine(double betaLine)
Returns an alpha line given a beta line.
Definition AlphaCube.h:109
double BetaSample(double alphaSample)
Returns a beta sample given an alpha sample.
Definition AlphaCube.h:145
double degrees() const
Get the angle in units of Degrees.
Definition Angle.h:232
@ Degrees
Degrees are generally considered more human readable, 0-360 is one circle, however most math does not...
Definition Angle.h:56
int Size() const
Accessor method that returns the number of bytes in the blob data.
Definition Blob.cpp:142
char * getBuffer()
Get the internal data buff of the Blob.
Definition Blob.cpp:546
PvlObject & Label()
Accessor method that returns a PvlObject containing the Blob label.
Definition Blob.cpp:151
QString getModelState() const
Get the CSM Model state string to re-create the CSM Model.
virtual void computeSolarLongitude(iTime et)
Computes the solar longitude for the given ephemeris time.
virtual double DetectorResolution()
Compute the detector resolution in meters per pixel for the current set point.
std::vector< int > getParameterIndices(csm::param::Set paramSet) const
Get the indices of the parameters that belong to a set.
CSMCamera(Cube &cube)
Constructor for an ISIS Camera model that uses a Community Sensor Model (CSM) for the principal trans...
Definition CSMCamera.cpp:62
void setTarget(Pvl label)
Set the Target object for the camera model.
QString getParameterUnits(int index)
Get the units of the parameter at a particular index.
void init(Cube &cube, QString pluginName, QString modelName, QString stateString)
Init method which performs most of the setup for the CSM Camera Model inside ISIS.
Definition CSMCamera.cpp:82
double getParameterCovariance(int index1, int index2)
Get the covariance between two parameters.
virtual bool SetImage(const double sample, const double line)
Set the image sample and line for the Camera Model and then compute the corresponding image time,...
virtual QList< QPointF > PixelIfovOffsets()
Returns the pixel ifov offsets from center of pixel.
virtual void setTime(const iTime &time)
Set the time and update the sensor position and orientation.
virtual double parentLine() const
Returns the currently set parent line for the camera model.
virtual SpicePosition * sunPosition() const
Get the SpicePosition object that contains the state information for the sun in J2000.
double getParameterValue(int index)
Get the value of a parameter.
virtual double SampleResolution()
Compute the sample resolution in meters per pixel for the current set point.
virtual void subSpacecraftPoint(double &lat, double &lon)
Get the latitude and longitude of the sub-spacecraft point at the currently set time.
virtual bool SetUniversalGround(const double latitude, const double longitude)
Set the latitude and longitude for the Camera Model and then compute the corresponding image time,...
virtual double LineResolution()
Compute the line resolution in meters per pixel for the current set point.
virtual double ObliqueLineResolution(bool useLocal=true)
Compute the oblique line resolution in meters per pixel for the current set point.
void csmToIsisPixel(csm::ImageCoord csmPixel, double &line, double &sample) const
Convert a CSM pixel coordinate to an ISIS pixel coordinate.
virtual double IncidenceAngle() const
Compute the incidence angle at the currently set ground point.
virtual double RightAscension()
Computes the Right Ascension of the currently set image coordinate.
csm::EcefCoord isisToCsmGround(const SurfacePoint &groundPt) const
Convert an ISIS ground point into a CSM ground point.
virtual std::vector< double > ImagePartials()
Compute the partial derivatives of the ground point with respect to the line and sample at the curren...
QString getParameterName(int index)
Get the name of the parameter.
virtual double CelestialNorthClockAngle()
Computes the celestial north clock angle at the current line/sample or ra/dec.
SurfacePoint csmToIsisGround(const csm::EcefCoord &groundPt) const
Convert a CSM ground point into an ISIS ground point.
virtual double Declination()
Computes the Declination of the currently set image coordinate.
void applyParameterCorrection(int index, double correction)
Adjust the value of a parameter.
virtual void subSolarPoint(double &lat, double &lon)
Returns the sub-solar latitude/longitude in universal coordinates (0-360 positive east,...
virtual std::vector< double > GroundPartials()
Compute the partial derivatives of the sample, line with respect to the x, y, z coordinates of the gr...
virtual bool SetRightAscensionDeclination(const double ra, const double dec)
Given the ra/dec compute the look direction.
virtual SpicePosition * instrumentPosition() const
Get the SpicePosition object the contains the state information for the sensor in J2000.
virtual double PhaseAngle() const
Compute the phase angle at the currently set ground point.
virtual double ObliqueSampleResolution(bool useLocal=true)
Compute the oblique sample resolution in meters per pixel for the current set point.
std::vector< double > sensorPositionBodyFixed() const
Get the position of the sensor in the body fixed coordinate system at the currently set time.
void isisToCsmPixel(double line, double sample, csm::ImageCoord &csmPixel) const
Convert an ISIS pixel coordinate to a CSM pixel coordinate.
virtual double EmissionAngle() const
Compute the emission angle at the currently set ground point.
virtual void instrumentBodyFixedPosition(double p[3]) const
Get the position of the sensor in the body fixed coordinate system at the currently set time.
virtual SpiceRotation * instrumentRotation() const
Get the SpiceRotation object the contains the orientation of the sensor relative to J2000.
virtual double parentSample() const
Returns the currently set parent sample for the camera model.
virtual double ObliqueDetectorResolution(bool useLocal=true)
Compute the oblique detector resolution in meters per pixel for the current set point.
virtual SpiceRotation * bodyRotation() const
Get the SpiceRotation object the contains the orientation of the target body relative to J2000.
virtual double SolarDistance() const
Computes the distance to the sun from the currently set ground point.
iTime m_refTime
CSM sensor model.
Definition CSMCamera.h:149
virtual double SlantDistance() const
Compute the slant distance from the sensor to the ground point at the currently set time.
virtual bool SetGround(Latitude latitude, Longitude longitude)
Set the latitude and longitude for the Camera Model and then compute the corresponding image time,...
virtual bool SetLookDirection(const std::vector< double > lookB)
Sets the look direction of the spacecraft.
virtual double targetCenterDistance() const
Calculates and returns the distance from the spacecraft to the target center at the currently set tim...
QString m_spacecraftNameLong
Full spacecraft name.
Definition Camera.h:499
AlphaCube * p_alphaCube
A pointer to the AlphaCube.
Definition Camera.h:504
virtual double Line() const
Returns the current line number.
Definition Camera.cpp:2740
double RaDecResolution()
Returns the RaDec resolution.
Definition Camera.cpp:1883
virtual double Sample() const
Returns the current sample number.
Definition Camera.cpp:2720
double p_childSample
Sample value for child.
Definition Camera.h:502
bool p_pointComputed
Flag showing if Sample/Line has been computed.
Definition Camera.h:506
QString m_instrumentNameShort
Shortened instrument name.
Definition Camera.h:498
QString m_spacecraftNameShort
Shortened spacecraft name.
Definition Camera.h:500
double p_childLine
Line value for child.
Definition Camera.h:503
QString m_instrumentNameLong
Full instrument name.
Definition Camera.h:497
IO Handler for Isis Cubes.
Definition Cube.h:168
void read(Blob &blob, const std::vector< PvlKeyword > keywords=std::vector< PvlKeyword >()) const
This method will read data from the specified Blob object.
Definition Cube.cpp:813
virtual QString fileName() const
Returns the opened cube's filename.
Definition Cube.cpp:1569
Pvl * label() const
Returns a pointer to the IsisLabel object associated with the cube.
Definition Cube.cpp:1707
Displacement is a signed length, usually in meters.
@ Kilometers
The distance is being specified in kilometers.
@ Meters
The distance is being specified in meters.
Distance measurement, usually in meters.
Definition Distance.h:34
bool isValid() const
Test if this distance has been initialized or not.
Definition Distance.cpp:192
@ Kilometers
The distance is being specified in kilometers.
Definition Distance.h:45
@ Meters
The distance is being specified in meters.
Definition Distance.h:43
Isis exception class.
Definition IException.h:91
@ User
A type of error that could only have occurred due to a mistake on the user's part (e....
Definition IException.h:126
@ Programmer
This error is for when a programmer made an API call that was illegal.
Definition IException.h:146
This class is designed to encapsulate the concept of a Latitude.
Definition Latitude.h:51
static Matrix pseudoinverse(const Matrix &matrix)
Returns the pseudoinverse of a matrix.
boost::numeric::ublas::matrix< double > Matrix
Definition for an Isis::LinearAlgebra::Matrix of doubles.
This class is designed to encapsulate the concept of a Longitude.
Definition Longitude.h:40
Container for cube-like labels.
Definition Pvl.h:119
Contains Pvl Groups and Pvl Objects.
Definition PvlObject.h:61
PvlKeyword & findKeyword(const QString &kname, FindOptions opts)
Finds a keyword in the current PvlObject, or deeper inside other PvlObjects and Pvlgroups within this...
virtual SurfacePoint GetSurfacePoint() const
Returns the surface point (most efficient accessor).
Definition Sensor.cpp:257
SpiceDouble m_lookB[3]
Look direction in body fixed.
Definition Sensor.h:234
bool HasSurfaceIntersection() const
Returns if the last call to either SetLookDirection or SetUniversalGround had a valid intersection wi...
Definition Sensor.cpp:188
bool m_newLookB
flag to indicate we need to recompute ra/dec
Definition Sensor.h:235
Distance LocalRadius() const
Returns the local radius at the intersection point.
Definition Sensor.cpp:269
Define shapes and provide utilities for Isis targets.
Definition ShapeModel.h:66
virtual void clearSurfacePoint()
Clears or resets the current surface point.
bool hasIntersection()
Returns intersection status.
void setHasIntersection(bool b)
Sets the flag to indicate whether this ShapeModel has an intersection.
QString name() const
Gets the shape name.
Target * m_target
Target of the observation.
Definition Spice.h:381
void radii(Distance r[3]) const
Returns the radii of the body in km.
Definition Spice.cpp:937
virtual Target * target() const
Returns a pointer to the target object.
Definition Spice.cpp:1380
iTime * m_et
Ephemeris time (read NAIF documentation for a detailed description)
Definition Spice.h:382
Obtain SPICE position information for a body.
Obtain SPICE rotation information for a body.
This class defines a body-fixed surface point.
Latitude GetLatitude() const
Return the body-fixed latitude for the surface point.
Longitude GetLongitude() const
Return the body-fixed longitude for the surface point.
This class is used to create and store valid Isis targets.
Definition Target.h:63
Parse and return pieces of a time string.
Definition iTime.h:65
This is free and unencumbered software released into the public domain.
Definition Apollo.h:16
const double DEG2RAD
Multiplier for converting from degrees to radians.
Definition Constants.h:43
const double Null
Value for an Isis Null pixel.
const double RAD2DEG
Multiplier for converting from radians to degrees.
Definition Constants.h:44
const double PI
The mathematical constant PI.
Definition Constants.h:40
Namespace for the standard library.