GLOCNavAlmDeltas.hpp
Go to the documentation of this file.
1 //==============================================================================
2 //
3 // This file is part of GNSSTk, the ARL:UT GNSS Toolkit.
4 //
5 // The GNSSTk is free software; you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published
7 // by the Free Software Foundation; either version 3.0 of the License, or
8 // any later version.
9 //
10 // The GNSSTk is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public
16 // License along with GNSSTk; if not, write to the Free Software Foundation,
17 // Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110, USA
18 //
19 // This software was developed by Applied Research Laboratories at the
20 // University of Texas at Austin.
21 // Copyright 2004-2022, The Board of Regents of The University of Texas System
22 //
23 //==============================================================================
24 
25 
26 //==============================================================================
27 //
28 // This software was developed by Applied Research Laboratories at the
29 // University of Texas at Austin, under contract to an agency or agencies
30 // within the U.S. Department of Defense. The U.S. Government retains all
31 // rights to use, duplicate, distribute, disclose, or release this software.
32 //
33 // Pursuant to DoD Directive 523024
34 //
35 // DISTRIBUTION STATEMENT A: This software has been approved for public
36 // release, distribution is unlimited.
37 //
38 //==============================================================================
39 #ifndef GNSSTK_GLOCNAVALMDELTAS_CPP
40 #define GNSSTK_GLOCNAVALMDELTAS_CPP
41 
45 #include "DebugTrace.hpp"
46 
47 namespace gnsstk
48 {
49  inline void GLOCNavAlm::Deltas ::
50  setData(double B, double Lk, const Uncorrected& uncor)
51  {
52  setdeltaa_a(B, Lk, uncor);
53  setdeltah(B, Lk, uncor);
54  setdeltal(B, Lk, uncor);
55  setdeltalambda(B, Lk, uncor);
56  setdeltai(B, Lk, uncor);
57  setdeltaLk(B, Lk, uncor);
58  }
59 
60 
61  inline void GLOCNavAlm::Deltas ::
62  setdeltaa_a(double B, double Lk, const Uncorrected& uncor)
63  {
65  2.0 * B * (1.0 - ((3.0/2.0)*sin(uncor.geti())*sin(uncor.geti())))*(uncor.getl()*cos(Lk)+uncor.geth()*sin(Lk))
66  + (B * sin(uncor.geti())*sin(uncor.geti())
67  * ((1.0/2.0)*uncor.geth()*sin(Lk)-((1.0/2.0)*uncor.getl()*cos(Lk))+cos(2.0*Lk)+
68  ((7.0/2.0)*uncor.getl()*cos(3.0*Lk))+((7.0/2.0)*uncor.geth()*sin(3.0*Lk))));
69  }
70 
71 
72  inline void GLOCNavAlm::Deltas ::
73  setdeltah(double B, double Lk, const Uncorrected& uncor)
74  {
75  Perturbations::h = B * (1.0 - ((3.0/2.0)*sin(uncor.geti())*sin(uncor.geti()))) *
76  (sin(Lk) + ((3.0/2.0)*uncor.getl()*sin(2.0*Lk)) - ((3.0/2.0)*uncor.geth()*cos(2.0*Lk)))
77  - ((1.0/4.0)*B*sin(uncor.geti())*sin(uncor.geti()) *
78  (sin(Lk)
79  -((7.0/3.0)*sin(3.0*Lk))
80  +(5.0*uncor.getl()*sin(2.0*Lk))
81  - ((17.0/2.0)*uncor.getl()*sin(4.0*Lk))
82  + ((17.0/2.0)*uncor.geth()*cos(4.0*Lk))
83  + uncor.geth()*cos(2.0*Lk)))
84  - ((1.0/2.0)*B*cos(uncor.geti())*cos(uncor.geti())*uncor.getl()*sin(2.0*Lk));
85  }
86 
87 
88  inline void GLOCNavAlm::Deltas ::
89  setdeltal(double B, double Lk, const Uncorrected& uncor)
90  {
91  Perturbations::l = B * (1.0 - ((3.0/2.0)*sin(uncor.geti())*sin(uncor.geti()))) *
92  (cos(Lk) + ((3.0/2.0)*uncor.getl()*cos(2.0*Lk)) + ((3.0/2.0)*uncor.geth()*sin(2.0*Lk)))
93  - ((1.0/4.0)*B*sin(uncor.geti())*sin(uncor.geti()) *
94  (-cos(Lk)
95  -((7.0/3.0)*cos(3.0*Lk))
96  -(5.0*uncor.geth()*sin(2.0*Lk))
97  -((17.0/2.0)*uncor.getl()*cos(4.0*Lk))
98  -((17.0/2.0)*uncor.geth()*sin(4.0*Lk))
99  +(uncor.getl()*cos(2.0*Lk))))
100  + ((1.0/2.0)*B*cos(uncor.geti())*cos(uncor.geti())*uncor.geth()*sin(2.0*Lk));
101  }
102 
103 
104  inline void GLOCNavAlm::Deltas ::
105  setdeltalambda(double B, double Lk, const Uncorrected& uncor)
106  {
107  Perturbations::lambda = -B * cos(uncor.geti()) *
108  (((7.0/2.0)*uncor.getl()*sin(Lk))
109  -((5.0/2.0)*uncor.geth()*cos(Lk))
110  -((1.0/2.0)*sin(2.0*Lk))
111  -((7.0/6.0)*uncor.getl()*sin(3.0*Lk))
112  +((7.0/6.0)*uncor.geth()*cos(3.0*Lk)));
113  }
114 
115 
116  inline void GLOCNavAlm::Deltas ::
117  setdeltai(double B, double Lk, const Uncorrected& uncor)
118  {
119  Perturbations::i = (1.0/2.0) * B * sin(uncor.geti()) * cos(uncor.geti()) *
120  ((-uncor.getl()*cos(Lk))
121  +(uncor.geth()*sin(Lk))
122  +(cos(2.0*Lk))
123  +((7.0/3.0)*uncor.getl()*cos(3.0*Lk))
124  +((7.0/3.0)*uncor.geth()*sin(3.0*Lk)));
125  }
126 
127 
128  inline void GLOCNavAlm::Deltas ::
129  setdeltaLk(double B, double Lk, const Uncorrected& uncor)
130  {
131  Perturbations::Lk = (2.0 * B * (1.0 - ((3.0/2.0)*sin(uncor.geti())*sin(uncor.geti()))) *
132  (((7.0/4.0)*uncor.getl()*sin(Lk)) - ((7.0/4.0)*uncor.geth()*cos(Lk))))
133  +((3.0*B*sin(uncor.geti())*sin(uncor.geti())) *
134  ((-(7.0/24.0)*uncor.geth()*cos(Lk))
135  -((7.0/24.0)*uncor.getl()*sin(Lk))
136  -((49.0/72.0)*uncor.geth()*cos(3.0*Lk))
137  +((49.0/72.0)*uncor.getl()*sin(3.0*Lk))
138  +((1.0/4.0)*sin(2.0*Lk))))
139  +((B*cos(uncor.geti())*cos(uncor.geti())) *
140  (((7.0/2.0)*uncor.getl()*sin(Lk))
141  -((5.0/2.0)*uncor.geth()*cos(Lk))
142  -((1.0/2.0)*sin(2.0*Lk))
143  -((7.0/6.0)*uncor.getl()*sin(3.0*Lk))
144  +((7.0/6.0)*uncor.geth()*cos(3.0*Lk))));
145  }
146 } // namespace gnsstk
147 
148 #endif //GNSSTK_GLOCNAVALMDELTAS_CPP
149 
gnsstk::GLOCNavAlm::Deltas::setdeltaLk
void setdeltaLk(double B, double Lk, const Uncorrected &uncor)
Definition: GLOCNavAlmDeltas.hpp:129
gnsstk::GLOCNavAlm::Perturbations::lambda
double lambda
longitude of ascending node
Definition: GLOCNavAlm.hpp:167
gnsstk::GLOCNavAlm::Perturbations::geti
double geti() const
Definition: GLOCNavAlm.hpp:161
gnsstk::GLOCNavAlm::Deltas::setdeltai
void setdeltai(double B, double Lk, const Uncorrected &uncor)
Definition: GLOCNavAlmDeltas.hpp:117
gnsstk::GLOCNavAlm::Deltas::setData
void setData(double B, double Lk, const Uncorrected &uncor)
Definition: GLOCNavAlmDeltas.hpp:50
gnsstk::GLOCNavAlm::Perturbations::l
double l
Quasi-Keplerian thing (true longitude?).
Definition: GLOCNavAlm.hpp:166
gnsstk
For Sinex::InputHistory.
Definition: BasicFramework.cpp:50
std::sin
double sin(gnsstk::Angle x)
Definition: Angle.hpp:144
gnsstk::GLOCNavAlm::Perturbations::Lk
double Lk
Mean longitude.
Definition: GLOCNavAlm.hpp:170
std::cos
double cos(gnsstk::Angle x)
Definition: Angle.hpp:146
gnsstk::GLOCNavAlm::Perturbations::a
double a
Semi-major axis.
Definition: GLOCNavAlm.hpp:164
gnsstk::GLOCNavAlm::Perturbations::geth
double geth() const
Definition: GLOCNavAlm.hpp:157
gnsstk::GLOCNavAlm::Deltas::setdeltal
void setdeltal(double B, double Lk, const Uncorrected &uncor)
Definition: GLOCNavAlmDeltas.hpp:89
gnsstk::GLOCNavAlm::Deltas::setdeltaa_a
void setdeltaa_a(double B, double Lk, const Uncorrected &uncor)
Definition: GLOCNavAlmDeltas.hpp:62
DebugTrace.hpp
gnsstk::GLOCNavAlm::Perturbations::getl
double getl() const
Definition: GLOCNavAlm.hpp:158
gnsstk::GLOCNavAlm::Perturbations::i
double i
inclination
Definition: GLOCNavAlm.hpp:169
gnsstk::GLOCNavAlm::Uncorrected
Definition: GLOCNavAlm.hpp:180
gnsstk::GLOCNavAlm::Perturbations::h
double h
Quasi-Keplerian thing (true latitude?).
Definition: GLOCNavAlm.hpp:165
gnsstk::GLOCNavAlm::Deltas::setdeltalambda
void setdeltalambda(double B, double Lk, const Uncorrected &uncor)
Definition: GLOCNavAlmDeltas.hpp:105
gnsstk::GLOCNavAlm::Deltas::setdeltah
void setdeltah(double B, double Lk, const Uncorrected &uncor)
Definition: GLOCNavAlmDeltas.hpp:73


gnsstk
Author(s):
autogenerated on Wed Oct 25 2023 02:40:39