File: observationPodAcceleration.cpp

package info (click to toggle)
groops 0%2Bgit20250907%2Bds-1
  • links: PTS, VCS
  • area: non-free
  • in suites: forky, sid
  • size: 11,140 kB
  • sloc: cpp: 135,607; fortran: 1,603; makefile: 20
file content (236 lines) | stat: -rw-r--r-- 9,787 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
/***********************************************/
/**
* @file observationPodAcceleration.cpp
*
* @brief Precise Orbit Data (POD) observations (Acceleration approach).
*
* @author Torsten Mayer-Guerr
* @date 2006-02-02
* update 2011-06-09 Norbert Zehentner
*
*/
/***********************************************/

#include "base/import.h"
#include "parallel/parallel.h"
#include "files/fileInstrument.h"
#include "files/fileSatelliteModel.h"
#include "misc/observation/observationMisc.h"
#include "classes/earthRotation/earthRotation.h"
#include "classes/parametrizationGravity/parametrizationGravity.h"
#include "classes/parametrizationAcceleration/parametrizationAcceleration.h"
#include "classes/observation/observation.h"
#include "classes/observation/observationPodAcceleration.h"

/***********************************************/

ObservationPodAcceleration::ObservationPodAcceleration(Config &config)
{
  try
  {
    FileName fileNameSatellite;
    FileName orbitName, starCameraName;
    UInt     epochCount, interpolationDegree;

    renameDeprecatedConfig(config, "satelliteModel", "inputfileSatelliteModel",     date2time(2020, 8, 19));
    renameDeprecatedConfig(config, "representation", "parametrizationGravity",      date2time(2020, 6,  3));
    renameDeprecatedConfig(config, "parameter",      "parametrizationAcceleration", date2time(2020, 6,  3));

    readConfig(config, "inputfileSatelliteModel",     fileNameSatellite,            Config::OPTIONAL, "{groopsDataDir}/satelliteModel/", "satellite macro model");
    readConfig(config, "rightHandSide",               rhs,                          Config::MUSTSET,  "",    "input for the reduced observation vector");
    readConfig(config, "inputfileOrbit",              orbitName,                    Config::MUSTSET,  "",    "used to evaluate the observation equations, not used as observations");
    readConfig(config, "inputfileStarCamera",         starCameraName,               Config::MUSTSET,  "",    "");
    readConfig(config, "earthRotation",               earthRotation,                Config::MUSTSET,  "",    "");
    readConfig(config, "ephemerides",                 ephemerides,                  Config::OPTIONAL, "jpl", "");
    readConfig(config, "parametrizationGravity",      parametrization,              Config::MUSTSET,  "",    "gravity field parametrization");
    readConfig(config, "parametrizationAcceleration", parametrizationAcceleration,  Config::DEFAULT,  "",    "orbit force parameters");
    readConfig(config, "interpolationDegree",         interpolationDegree,          Config::DEFAULT,  "8",   "orbit differentation  by polynomial approximation of degree n");
    readConfig(config, "numberOfEpochs",              epochCount,                   Config::DEFAULT,  "9",   "number of used Epochs for polynom computation");
    readConfig(config, "covariancePod",               covPod,                       Config::OPTIONAL, "",    "covariance matrix of kinematic orbits");
    if(isCreateSchema(config)) return;

    // test if epochCount and polynomial degree are correct
    // ---------------------
    if(interpolationDegree%2 != 0)
      throw(Exception("polnomial degree for interpolation must be even."));
    if(epochCount <= interpolationDegree)
      throw(Exception("Epochs must be greater than interpolationdegree, at least by one."));
    if(epochCount%2 != 1)
      throw(Exception("Epochs must be odd (symmetry around central point)."));

    if(!fileNameSatellite.empty())
      readFileSatelliteModel(fileNameSatellite, satellite);

    // test instrument files
    // ---------------------
    orbitFile.open(orbitName);
    starCameraFile.open(starCameraName);
    InstrumentFile::checkArcCount({orbitFile, starCameraFile});
    for(UInt j=0; j<rhs.size(); j++)
      InstrumentFile::checkArcCount({orbitFile, *rhs.at(j)->orbitFile, *rhs.at(j)->accelerometerFile});

    // orbit differentation coefficients (position -> acceleration) by means of QR-decomposition
    // with or without overdetermination
    // -----------------------
    Matrix P(epochCount, interpolationDegree+1);
    for(UInt i=0; i<epochCount; i++)
      for(UInt n=0; n<=interpolationDegree; n++)
        P(i,n) = ((n==0) ? 1.0 : std::pow(i-(epochCount-1.)/2., n));

    const Vector tau = QR_decomposition(P);
    coeff = Vector(epochCount);
    coeff(2) = 2.0;
    triangularSolve(1., P.row(0,interpolationDegree+1).trans(), coeff);
    QMult(P, tau, coeff);
  }
  catch(std::exception &e)
  {
    GROOPS_RETHROW(e)
  }
}

/***********************************************/

Bool ObservationPodAcceleration::setInterval(const Time &timeStart, const Time &timeEnd)
{
  try
  {
    Bool change = FALSE;
    change = parametrization->setInterval(timeStart, timeEnd)             || change;
    change = parametrizationAcceleration->setInterval(timeStart, timeEnd) || change;
    return change;
  }
  catch(std::exception &e)
  {
    GROOPS_RETHROW(e)
  }
}

/***********************************************/

UInt ObservationPodAcceleration::parameterCount() const
{
  return parametrization->parameterCount() + parametrizationAcceleration->parameterCount();
}

/***********************************************/

void ObservationPodAcceleration::parameterName(std::vector<ParameterName> &name) const
{
  if(parametrization)
    parametrization->parameterName(name);

  if(parametrizationAcceleration)
  {
    parametrizationAcceleration->parameterName(name);
    const std::string satelliteName = satellite ? satellite->satelliteName : "satellite";
    for(UInt i=name.size(); i-->name.size()-parametrizationAcceleration->parameterCount();)
      name.at(i).object = satelliteName;
  }
}

/***********************************************/

void ObservationPodAcceleration::observation(UInt arcNo, Matrix &l, Matrix &A, Matrix &B)
{
  try
  {
    OrbitArc      orbit      = orbitFile.readArc(arcNo);
    StarCameraArc starCamera = starCameraFile.readArc(arcNo);
    UInt          rhsCount   = rhs.size();
    UInt          posCount   = orbit.size();
    UInt          obsCount   = posCount-coeff.rows()+1;
    UInt          half       = (coeff.rows()-1)/2;
    Double        dt         = medianSampling(orbit.times()).seconds();
    Arc::checkSynchronized({orbit, starCamera});

    parametrizationAcceleration->setIntervalArc(orbit.at(half).time, orbit.at(obsCount+half-1).time+seconds2time(dt));

    // calculate earthrotation
    // -----------------------
    std::vector<Rotary3d> rotEarth(obsCount);
    for(UInt k=0; k<obsCount; k++)
      rotEarth.at(k) = earthRotation->rotaryMatrix(orbit.at(k+half).time);

    // reference acceleration
    // ----------------------
    std::vector<std::vector<Vector3d>> reference(rhsCount);
    for(UInt j=0; j<rhsCount; j++)
    {
      reference.at(j).resize(obsCount);
      AccelerometerArc accelerometer = rhs.at(j)->accelerometerFile->readArc(arcNo);
      Arc::checkSynchronized({orbit, accelerometer});

      for(UInt k=0; k<obsCount; k++)
      {
        Vector3d g = rhs.at(j)->forces->acceleration(satellite, orbit.at(k+half).time, orbit.at(k+half).position, orbit.at(k+half).velocity,
                                                     starCamera.at(k+half).rotary, rotEarth.at(k), earthRotation, ephemerides);
        if(accelerometer.size())
          g += rotEarth.at(k).rotate(starCamera.at(k+half).rotary.rotate(accelerometer.at(k+half).acceleration)); // accelerometer
        reference.at(j).at(k) = g; // reference observations
      }
    }

    // reduced observations
    // --------------------
    l = Matrix(3*obsCount, rhsCount);
    for(UInt j=0; j<rhsCount; j++)
    {
      OrbitArc orbitPod = rhs.at(j)->orbitFile->readArc(arcNo);
      for(UInt k=0; k<obsCount; k++)
        {
          Vector3d a;
          for(UInt i=0; i<coeff.rows(); i++)
            a += (coeff(i)/dt/dt) * orbitPod.at(k+i).position;
          a = rotEarth.at(k).rotate(a);

          l(3*k+0,j) = a.x() - reference.at(j).at(k).x();
          l(3*k+1,j) = a.y() - reference.at(j).at(k).y();
          l(3*k+2,j) = a.z() - reference.at(j).at(k).z();
        }
    }

    // Design matrix A and B
    // ---------------------
    A = Matrix(3*obsCount, parametrization->parameterCount() + parametrizationAcceleration->parameterCount());
    B = Matrix(3*obsCount, parametrizationAcceleration->parameterCountArc());
    // gravity
    for(UInt k=0; k<obsCount; k++)
      parametrization->gravity(orbit.at(k+half).time, rotEarth.at(k).rotate(orbit.at(k+half).position), A.slice(3*k, 0, 3, parametrization->parameterCount()));

    // orbit parameters
    for(UInt k=0; k<obsCount; k++)
    {
      parametrizationAcceleration->compute(satellite, orbit.at(k+half).time, orbit.at(k+half).position, orbit.at(k+half).velocity,
                                           starCamera.at(k+half).rotary, rotEarth.at(k), ephemerides,
                                           A.slice(3*k, parametrization->parameterCount(), 3, parametrizationAcceleration->parameterCount()), B.row(3*k,3));
    }

    // decorrelation
    // -------------
    if(covPod)
    {
      // linearized variance propagation
      Matrix D(3*obsCount, 3*posCount);
      for(UInt k=0; k<obsCount; k++)
        for(UInt i=0; i<coeff.rows(); i++)
          axpy((coeff(i)/dt/dt), rotEarth.at(k).matrix(), D.slice(3*k,3*(k+i),3,3));
      Matrix C = covPod->covariance(arcNo, orbit);
      Matrix DCD = D * C * D.trans();
      DCD.setType(Matrix::SYMMETRIC);
      cholesky(DCD);

      // apply Cholesky matrix
      triangularSolve(1., DCD.trans(), A);
      triangularSolve(1., DCD.trans(), l);
      if(B.size())
        triangularSolve(1., DCD.trans(), B);
    }
  }
  catch(std::exception &e)
  {
    GROOPS_RETHROW(e)
  }
}

/***********************************************/