-
Notifications
You must be signed in to change notification settings - Fork 73
Чебышёв-Picard Iteration #4409
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
base: master
Are you sure you want to change the base?
Чебышёв-Picard Iteration #4409
Changes from all commits
dfecc88
b95ba19
d3585e7
31bb19c
6143f45
7103c72
e908836
5e1cf93
f03fd4f
5f9bb7b
2bc7369
c4429b7
72b52df
cb02b0e
37e426b
7895c44
43eef1b
3c22631
a73c0d2
0192266
b957136
ac29708
2849aa2
6c5925a
afe326d
8c9c620
1965e8d
a346203
1a66062
6bbbaf2
16bcf39
56a2a2e
903079c
a11e5b5
0f235d7
fed4c5f
ec46317
17be53f
bfb2881
d320f0e
0139a50
abf6275
697bf9a
6b30b63
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
|
|
@@ -44,6 +44,28 @@ @article{AriasCharlotFeisselLestrade1995 | |
| volume = {303}, | ||
| } | ||
|
|
||
| @article{BaiJunkins2010, | ||
| author = {Bai, Xiaoli and Junkins, John L}, | ||
| date = {2010}, | ||
| journaltitle = {Advances in the Astronautical Sciences}, | ||
| number = {2}, | ||
| pages = {1459--1476}, | ||
| title = {Solving initial value problems by the Picard-Chebyshev method with NVIDIA GPUs}, | ||
| volume = {136}, | ||
| } | ||
|
|
||
| @article{BaiJunkins2012, | ||
| author = {Bai, Xiaoli and Junkins, John L}, | ||
| publisher = {Springer}, | ||
| date = {2012}, | ||
| doi = {10.1007/s40295-013-0021-6}, | ||
| journaltitle = {The Journal of the Astronautical Sciences}, | ||
| number = {1}, | ||
| pages = {327--351}, | ||
| title = {Modified Chebyshev-Picard iteration methods for solution of initial value problems}, | ||
| volume = {59}, | ||
| } | ||
|
|
||
| @article{Bateman1938, | ||
| author = {Bateman, Harry}, | ||
| date = {1938-01}, | ||
|
|
@@ -393,6 +415,18 @@ @article{ChinKidwell2000 | |
| volume = {62}, | ||
| } | ||
|
|
||
| @article{ClenshawNorton1963, | ||
| author = {Clenshaw, C. W. and Norton, H. J.}, | ||
| publisher = {Oxford University Press}, | ||
| date = {1963}, | ||
| doi = {10.1093/comjnl/6.1.88}, | ||
| journaltitle = {The Computer Journal}, | ||
| number = {1}, | ||
| pages = {88--92}, | ||
| title = {The solution of nonlinear ordinary differential equations in Chebyshev series}, | ||
| volume = {6}, | ||
| } | ||
|
|
||
| @article{Cohen1951, | ||
| author = {Cohen, Jr., A. C.}, | ||
| date = {1951-06}, | ||
|
|
@@ -581,6 +615,18 @@ @article{FarrésLaskarBlanesCasasMakazagaMurua2013 | |
| volume = {116}, | ||
| } | ||
|
|
||
| @article{FeaginNacozy1983, | ||
| author = {Feagin, Terry and Nacozy, Paul}, | ||
| publisher = {Springer}, | ||
| date = {1983}, | ||
| doi = {10.1007/BF01232802}, | ||
| journaltitle = {Celestial mechanics}, | ||
| number = {2}, | ||
| pages = {107--115}, | ||
| title = {Matrix formulation of the Picard method for parallel computation}, | ||
| volume = {29}, | ||
| } | ||
|
|
||
| @article{FiengaMacheLaskarGastineau2008, | ||
| author = {Fienga, A. and Manche, H. and Laskar, J. and Gastineau, M.}, | ||
| url = {http://dx.doi.org/10.1051/0004-6361:20066607}, | ||
|
|
@@ -654,6 +700,13 @@ @article{Fornberg1988 | |
| volume = {51}, | ||
| } | ||
|
|
||
| @article{Fukushima1997, | ||
| author = {Fukushima, Toshio}, | ||
| doi = {10.1086/118404}, | ||
| journaltitle = {Astronomical Journal v. 113, p. 1909-1914 (1997)}, | ||
rnlahaye marked this conversation as resolved.
Show resolved
Hide resolved
|
||
| title = {Picard iteration method, Chebyshev polynomial approximation, and global numerical integration of dynamical motions}, | ||
| } | ||
rnlahaye marked this conversation as resolved.
Show resolved
Hide resolved
|
||
|
|
||
| @article{Fukushima2009a, | ||
| author = {Fukushima, Toshio}, | ||
| date = {2009-10-25}, | ||
|
|
@@ -2242,6 +2295,14 @@ @report{BlancSchgounn1996 | |
| title = {AVISO/Altimetry, AVISO User Handbook for Merged TOPEX/\-POSEIDON products, Edition 3.0}, | ||
| } | ||
|
|
||
| @report{ByrdLee1969, | ||
| author = {Byrd, P. F. and Lee, K. L.}, | ||
| date = {1969}, | ||
| number = {19690027545}, | ||
| title = {Chebyshev series solution of nonlinear ordinary differential equations-Initial-value problems}, | ||
| type = {techreport}, | ||
rnlahaye marked this conversation as resolved.
Show resolved
Hide resolved
|
||
| } | ||
|
|
||
| @report{HairerWanner2010, | ||
| author = {Hairer, Ernst and Wanner, Gerhard}, | ||
| institution = {Department of Mathematics, University of Geneva, Switzerland}, | ||
|
|
@@ -2393,6 +2454,14 @@ @thesis{Chandler73 | |
| type = {phdthesis}, | ||
| } | ||
|
|
||
| @thesis{Macomber2015, | ||
| author = {Macomber, Brent David}, | ||
| url = {https://oaktrust.library.tamu.edu/server/api/core/bitstreams/e94bc2a8-7916-4d02-9aef-e7493a07ccb9/content}, | ||
| date = {2015}, | ||
| title = {Enhancements to Chebyshev-Picard iteration efficiency for generally perturbed orbits and constrained dynamical systems}, | ||
| type = {phdthesis}, | ||
| } | ||
|
Member
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Same, this was eventually published, see https://www.researchgate.net/publication/320007623_Enhancements_to_modified_Chebyshev-Picard_iteration_efficiency_for_perturbed_orbit_propagation.
Contributor
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. While the published paper contains the second-order vector/matrix formulation, it does not contain the first-order one (which is what I'm using here). I added a URL to the citation. |
||
|
|
||
| @unpublished{Andrew2017, | ||
| author = {“Andrew”}, | ||
| date = {2017-05-10}, | ||
|
|
||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,149 @@ | ||
| // The files containing the tree of of child classes of `Integrator` must be | ||
| // included in the order of inheritance to avoid circular dependencies. This | ||
| // class will end up being reincluded as part of the implementation of its | ||
| // parent. | ||
| #ifndef PRINCIPIA_INTEGRATORS_INTEGRATORS_HPP_ | ||
| #include "integrators/integrators.hpp" | ||
| #else | ||
| #ifndef PRINCIPIA_INTEGRATORS_ЧЕБЫШЁВ_PICARD_ITERATOR_HPP_ | ||
| #define PRINCIPIA_INTEGRATORS_ЧЕБЫШЁВ_PICARD_ITERATOR_HPP_ | ||
|
|
||
| #include <memory> | ||
| #include <vector> | ||
|
|
||
| #include "absl/status/status.h" | ||
| #include "base/not_null.hpp" | ||
| #include "integrators/ordinary_differential_equations.hpp" | ||
| #include "numerics/fixed_arrays.hpp" | ||
| #include "quantities/quantities.hpp" | ||
| #include "serialization/integrators.pb.h" | ||
|
|
||
| namespace principia { | ||
| namespace integrators { | ||
| namespace _чебышёв_picard_iterator { | ||
| namespace internal { | ||
|
|
||
| using namespace principia::base::_not_null; | ||
| using namespace principia::integrators::_integrators; | ||
|
Member
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Let's have an include for this.
Contributor
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. The include is near the top, line 6 (like the other integrators). |
||
| using namespace principia::integrators::_ordinary_differential_equations; | ||
| using namespace principia::numerics::_fixed_arrays; | ||
| using namespace principia::quantities::_quantities; | ||
|
|
||
| struct ЧебышёвPicardIterationParams { | ||
|
Member
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. So we generally distinguish the method (see For instance for the ephemeris, the method is It's clear that
Contributor
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. The working matrices are fixed once I see that all the methods in (in alternatively one could have (or am I misunderstanding?)
Member
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I think
Contributor
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Ok, added This enabled me to switch from UnboundedMatrix to FixedMatrix, which resulted in me finding and fixing a bug in the M > N case. |
||
| // The maximum allowed number of Picard iterations per step. If iteration has | ||
| // not stopped (according to the stopping criterion) by the final step, the | ||
| // iteration will be considered to have diverged. | ||
| std::int64_t max_iterations; | ||
|
|
||
| // If the maximum absolute difference between successive state approximations | ||
| // is less than this for two Picard iterations in a row, iteration will be | ||
| // considered to have converged. | ||
| double stopping_criterion; | ||
| }; | ||
|
|
||
| // This class solves ordinary differential equations of the form x′ = f(x, t) | ||
| // using Чебышёв-Picard iteration. | ||
| // | ||
| // Чебышёв-Picard iteration combines Чебышёв interpolation with Picard | ||
| // iteration; it was first proposed in [CN63]. It works as follows: | ||
| // | ||
| // * Start with some initial approximation x⁰ (the constant function | ||
| // x⁰(t) = x(t₀) is typically used). | ||
| // * Repeatedly construct xⁱ⁺¹ from xⁱ by | ||
| // * Approximating f(t, xⁱ) using Чебышёв interpolation. | ||
| // * Approximating x₀ + ∫ₜ₀ᵗ f(t, xⁱ(t)) dt by integrating the above | ||
| // approximation. This will be xⁱ⁺¹. | ||
| // * The sequence xⁱ should converge to (a Чебышёв approximation of) x within | ||
| // some interval of t₀. | ||
| // | ||
| // Due to the properties of Чебышёв polynomials, the iteration may be peformed | ||
| // simply using linear algebra operations (in addition to the required | ||
| // evaluation of f at Чебышёв nodes) [FN83]. | ||
| // | ||
| // This code uses the formulation from [Mac15]. | ||
|
|
||
| template<typename Method, typename ODE_> | ||
| class ЧебышёвPicardIterator : public FixedStepSizeIntegrator<ODE_> { | ||
|
Member
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I am a bit unhappy that this is called (If we decide to change the name, let's do it last.)
Contributor
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Happy to change the name after everything else is resolved. |
||
| public: | ||
| using ODE = ODE_; | ||
| using AppendState = typename Integrator<ODE>::AppendState; | ||
|
|
||
| class Instance : public FixedStepSizeIntegrator<ODE>::Instance { | ||
| public: | ||
| absl::Status Solve(ODE::IndependentVariable const& t_final) override; | ||
|
|
||
| ЧебышёвPicardIterator const& integrator() const override; | ||
|
|
||
| not_null<std::unique_ptr<typename Integrator<ODE>::Instance>> Clone() | ||
| const override; | ||
|
|
||
| private: | ||
| Instance(InitialValueProblem<ODE> const& problem, | ||
rnlahaye marked this conversation as resolved.
Show resolved
Hide resolved
|
||
| AppendState const& append_state, | ||
| Time const& step, | ||
| ЧебышёвPicardIterator const& integrator); | ||
|
|
||
| // The dimension of the ODE. | ||
| static constexpr std::int64_t n = | ||
| std::tuple_size_v<typename ODE::DependentVariables>; | ||
|
|
||
| ЧебышёвPicardIterator const& integrator_; | ||
|
|
||
| // Working variables which are stored here so we don't need to reallocate | ||
| // them on each Solve call. | ||
|
|
||
| // Stores the nodes rescaled to the current step. | ||
| std::vector<typename ODE::IndependentVariable> t_; | ||
|
|
||
| // Controls the boundary condition. | ||
| FixedMatrix<double, Method::M + 1, n> CₓX₀_; | ||
|
|
||
| // Xⁱ is an (M + 1)×n matrix containing the values of the dependent | ||
| // variables at each node. | ||
| FixedMatrix<double, Method::M + 1, n> Xⁱ_; | ||
| FixedMatrix<double, Method::M + 1, n> Xⁱ⁺¹_; | ||
|
|
||
| // The computed derivative (at each node, for the current iteration). | ||
| FixedMatrix<double, Method::M + 1, n> yʹ_; | ||
|
|
||
| friend class ЧебышёвPicardIterator; | ||
| }; | ||
|
|
||
| // Constructs a ЧебышёвPicardIterator with the given parameters. | ||
| explicit ЧебышёвPicardIterator(ЧебышёвPicardIterationParams const& params); | ||
|
Member
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. As mentioned above, some of the things in
Contributor
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. How does this work with |
||
|
|
||
| ЧебышёвPicardIterationParams const& params() const; | ||
|
|
||
| not_null<std::unique_ptr<typename Integrator<ODE>::Instance>> NewInstance( | ||
| InitialValueProblem<ODE> const& problem, | ||
| AppendState const& append_state, | ||
| Time const& step) const override; | ||
|
|
||
| void WriteToMessage( | ||
| not_null<serialization::FixedStepSizeIntegrator*> message) const override; | ||
|
|
||
| private: | ||
| ЧебышёвPicardIterationParams params_; | ||
|
|
||
| // The nodes used for function evaluation. | ||
| // | ||
| // These are Чебышёв nodes of the second kind. | ||
| FixedVector<double, Method::M + 1> nodes_; | ||
|
|
||
| // The product of 1.31a and 1.31b from [Mac15]. | ||
| FixedMatrix<double, Method::M + 1, Method::M + 1> CₓCα_; | ||
| }; | ||
|
|
||
| } // namespace internal | ||
|
|
||
| using internal::ЧебышёвPicardIterationParams; | ||
| using internal::ЧебышёвPicardIterator; | ||
|
|
||
| } // namespace _чебышёв_picard_iterator | ||
| } // namespace integrators | ||
| } // namespace principia | ||
|
|
||
| #include "integrators/чебышёв_picard_iterator_body.hpp" | ||
|
|
||
| #endif // PRINCIPIA_INTEGRATORS_ЧЕБЫШЁВ_PICARD_ITERATOR_HPP_ | ||
| #endif // PRINCIPIA_INTEGRATORS_INTEGRATORS_HPP_ | ||
Uh oh!
There was an error while loading. Please reload this page.