not_null<std::unique_ptr<Ephemeris<Frame>>> Ephemeris<Frame>::ReadFromMessage( serialization::Ephemeris const& message) { std::vector<not_null<std::unique_ptr<MassiveBody const>>> bodies; for (auto const& body : message.body()) { bodies.push_back(MassiveBody::ReadFromMessage(body)); } auto const fitting_tolerance = Length::ReadFromMessage(message.fitting_tolerance()); std::unique_ptr<FixedStepParameters> parameters; bool const is_pre_буняковский = message.has_planetary_integrator(); if (is_pre_буняковский) { auto const& planetary_integrator = FixedStepSizeIntegrator<NewtonianMotionEquation>::ReadFromMessage( message.planetary_integrator()); CHECK(message.has_step()); auto const step = Time::ReadFromMessage(message.step()); parameters = std::make_unique<FixedStepParameters>(planetary_integrator, step); } else { parameters.reset(new FixedStepParameters( FixedStepParameters::ReadFromMessage(message.fixed_step_parameters()))); } // Dummy initial state and time. We'll overwrite them later. std::vector<DegreesOfFreedom<Frame>> const initial_state( bodies.size(), DegreesOfFreedom<Frame>(Position<Frame>(), Velocity<Frame>())); Instant const initial_time; auto ephemeris = make_not_null_unique<Ephemeris<Frame>>( std::move(bodies), initial_state, initial_time, fitting_tolerance, *parameters); ephemeris->last_state_ = NewtonianMotionEquation::SystemState::ReadFromMessage( message.last_state()); int index = 0; ephemeris->bodies_to_trajectories_.clear(); ephemeris->trajectories_.clear(); for (auto const& trajectory : message.trajectory()) { not_null<MassiveBody const*> const body = ephemeris->bodies_[index].get(); not_null<std::unique_ptr<ContinuousTrajectory<Frame>>> deserialized_trajectory = ContinuousTrajectory<Frame>::ReadFromMessage(trajectory); ephemeris->trajectories_.push_back(deserialized_trajectory.get()); ephemeris->bodies_to_trajectories_.emplace( body, std::move(deserialized_trajectory)); ++index; } if (message.has_t_max()) { ephemeris->checkpoints_.push_back(ephemeris->GetCheckpoint()); ephemeris->Prolong(Instant::ReadFromMessage(message.t_max())); } return ephemeris; }
void Ephemeris<Frame>::FlowWithFixedStep( std::vector<not_null<DiscreteTrajectory<Frame>*>> const& trajectories, std::vector<IntrinsicAcceleration> const& intrinsic_accelerations, Instant const& t, FixedStepParameters const& parameters) { VLOG(1) << __FUNCTION__ << " " << NAMED(parameters.step_) << " " << NAMED(t); if (empty() || t > t_max()) { Prolong(t); } std::vector<typename ContinuousTrajectory<Frame>::Hint> hints(bodies_.size()); NewtonianMotionEquation massless_body_equation; massless_body_equation.compute_acceleration = std::bind(&Ephemeris::ComputeMasslessBodiesTotalAccelerations, this, std::cref(intrinsic_accelerations), _1, _2, _3, &hints); typename NewtonianMotionEquation::SystemState initial_state; for (auto const& trajectory : trajectories) { auto const trajectory_last = trajectory->last(); auto const last_degrees_of_freedom = trajectory_last.degrees_of_freedom(); // TODO(phl): why do we keep rewriting this? Should we check consistency? initial_state.time = trajectory_last.time(); initial_state.positions.push_back(last_degrees_of_freedom.position()); initial_state.velocities.push_back(last_degrees_of_freedom.velocity()); } IntegrationProblem<NewtonianMotionEquation> problem; problem.equation = massless_body_equation; #if defined(WE_LOVE_228) typename NewtonianMotionEquation::SystemState last_state; problem.append_state = [&last_state]( typename NewtonianMotionEquation::SystemState const& state) { last_state = state; }; #else problem.append_state = std::bind(&Ephemeris::AppendMasslessBodiesState, _1, std::cref(trajectories)); #endif problem.t_final = t; problem.initial_state = &initial_state; parameters.integrator_->Solve(problem, parameters.step_); #if defined(WE_LOVE_228) // The |positions| are empty if and only if |append_state| was never called; // in that case there was not enough room to advance the |trajectories|. if (!last_state.positions.empty()) { AppendMasslessBodiesState(last_state, trajectories); } #endif }
void BM_BarycentricRotatingDynamicFrame( benchmark::State& state) { // NOLINT(runtime/references) Time const Δt = 5 * Minute; int const steps = state.range_x(); SolarSystem<ICRFJ2000Equator> solar_system; solar_system.Initialize( SOLUTION_DIR / "astronomy" / "gravity_model.proto.txt", SOLUTION_DIR / "astronomy" / "initial_state_jd_2433282_500000000.proto.txt"); auto const ephemeris = solar_system.MakeEphemeris( McLachlanAtela1992Order5Optimal<Position<ICRFJ2000Equator>>(), 45 * Minute, 5 * Milli(Metre)); ephemeris->Prolong(solar_system.epoch() + steps * Δt); not_null<MassiveBody const*> const earth = solar_system.massive_body(*ephemeris, "Earth"); not_null<MassiveBody const*> const venus = solar_system.massive_body(*ephemeris, "Venus"); MasslessBody probe; Position<ICRFJ2000Equator> probe_initial_position = ICRFJ2000Equator::origin + Displacement<ICRFJ2000Equator>( {0.5 * AstronomicalUnit, -1 * AstronomicalUnit, 0 * AstronomicalUnit}); Velocity<ICRFJ2000Equator> probe_velocity = Velocity<ICRFJ2000Equator>({0 * SIUnit<Speed>(), 100 * Kilo(Metre) / Second, 0 * SIUnit<Speed>()}); DiscreteTrajectory<ICRFJ2000Equator> probe_trajectory; FillLinearTrajectory<ICRFJ2000Equator, DiscreteTrajectory>( probe_initial_position, probe_velocity, solar_system.epoch(), Δt, steps, &probe_trajectory); BarycentricRotatingDynamicFrame<ICRFJ2000Equator, Rendering> dynamic_frame(ephemeris.get(), earth, venus); while (state.KeepRunning()) { auto v = ApplyDynamicFrame(&probe, &dynamic_frame, probe_trajectory.Begin(), probe_trajectory.End()); } }
std::unique_ptr<Ephemeris<Frame>> Ephemeris<Frame>::ReadFromPreBourbakiMessages( google::protobuf::RepeatedPtrField< serialization::Plugin::CelestialAndProperties> const& messages, Length const& fitting_tolerance, typename Ephemeris<Frame>::FixedStepParameters const& fixed_parameters) { LOG(INFO) << "Reading "<< messages.SpaceUsedExcludingSelf() << " bytes in pre-Bourbaki compatibility mode "; std::vector<not_null<std::unique_ptr<MassiveBody const>>> bodies; std::vector<DegreesOfFreedom<Frame>> initial_state; std::vector<std::unique_ptr<DiscreteTrajectory<Frame>>> histories; std::set<Instant> initial_time; std::set<Instant> final_time; for (auto const& message : messages) { serialization::Celestial const& celestial = message.celestial(); bodies.emplace_back(MassiveBody::ReadFromMessage(celestial.body())); histories.emplace_back(DiscreteTrajectory<Frame>::ReadFromMessage( celestial.history_and_prolongation().history(), /*forks=*/{})); auto const prolongation = DiscreteTrajectory<Frame>::ReadPointerFromMessage( celestial.history_and_prolongation().prolongation(), histories.back().get()); typename DiscreteTrajectory<Frame>::Iterator const history_begin = histories.back()->Begin(); initial_state.push_back(history_begin.degrees_of_freedom()); initial_time.insert(history_begin.time()); final_time.insert(prolongation->last().time()); } CHECK_EQ(1, initial_time.size()); CHECK_EQ(1, final_time.size()); LOG(INFO) << "Initial time is " << *initial_time.cbegin() << ", final time is " << *final_time.cbegin(); // Construct a new ephemeris using the bodies and initial states and time // extracted from the serialized celestials. auto ephemeris = std::make_unique<Ephemeris<Frame>>(std::move(bodies), initial_state, *initial_time.cbegin(), fitting_tolerance, fixed_parameters); // Extend the continuous trajectories using the data from the discrete // trajectories. std::set<Instant> last_state_time; for (int i = 0; i < histories.size(); ++i) { not_null<MassiveBody const*> const body = ephemeris->unowned_bodies_[i]; auto const& history = histories[i]; int const j = ephemeris->serialization_index_for_body(body); auto continuous_trajectory = ephemeris->trajectories_[j]; typename DiscreteTrajectory<Frame>::Iterator it = history->Begin(); Instant last_time = it.time(); DegreesOfFreedom<Frame> last_degrees_of_freedom = it.degrees_of_freedom(); for (; it != history->End(); ++it) { Time const duration_since_last_time = it.time() - last_time; if (duration_since_last_time == fixed_parameters.step_) { // A time in the discrete trajectory that is aligned on the continuous // trajectory. last_time = it.time(); last_degrees_of_freedom = it.degrees_of_freedom(); continuous_trajectory->Append(last_time, last_degrees_of_freedom); } else if (duration_since_last_time > fixed_parameters.step_) { // A time in the discrete trajectory that is not aligned on the // continuous trajectory. Stop here, we'll use prolong to recompute the // rest. break; } } // Fill the |last_state_| for this body. It will be the starting state for // Prolong. last_state_time.insert(last_time); ephemeris->last_state_.positions[j] = last_degrees_of_freedom.position(); ephemeris->last_state_.velocities[j] = last_degrees_of_freedom.velocity(); } CHECK_EQ(1, last_state_time.size()); ephemeris->last_state_.time = *last_state_time.cbegin(); LOG(INFO) << "Last time in discrete trajectories is " << *last_state_time.cbegin(); // Prolong the ephemeris to the final time. This might create discrepancies // from the discrete trajectories. ephemeris->Prolong(*final_time.cbegin()); return ephemeris; }
bool Ephemeris<Frame>::FlowWithAdaptiveStep( not_null<DiscreteTrajectory<Frame>*> const trajectory, IntrinsicAcceleration intrinsic_acceleration, Instant const& t, AdaptiveStepParameters const& parameters, std::int64_t const max_ephemeris_steps) { Instant const& trajectory_last_time = trajectory->last().time(); if (trajectory_last_time == t) { return true; } std::vector<not_null<DiscreteTrajectory<Frame>*>> const trajectories = {trajectory}; std::vector<IntrinsicAcceleration> const intrinsic_accelerations = {std::move(intrinsic_acceleration)}; // The |min| is here to prevent us from spending too much time computing the // ephemeris. The |max| is here to ensure that we always try to integrate // forward. We use |last_state_.time.value| because this is always finite, // contrary to |t_max()|, which is -∞ when |empty()|. Instant const t_final = std::min(std::max(last_state_.time.value + max_ephemeris_steps * parameters_.step(), trajectory_last_time + parameters_.step()), t); Prolong(t_final); std::vector<typename ContinuousTrajectory<Frame>::Hint> hints(bodies_.size()); NewtonianMotionEquation massless_body_equation; massless_body_equation.compute_acceleration = std::bind(&Ephemeris::ComputeMasslessBodiesTotalAccelerations, this, std::cref(intrinsic_accelerations), _1, _2, _3, &hints); typename NewtonianMotionEquation::SystemState initial_state; auto const trajectory_last = trajectory->last(); auto const last_degrees_of_freedom = trajectory_last.degrees_of_freedom(); initial_state.time = trajectory_last.time(); initial_state.positions.push_back(last_degrees_of_freedom.position()); initial_state.velocities.push_back(last_degrees_of_freedom.velocity()); IntegrationProblem<NewtonianMotionEquation> problem; problem.equation = massless_body_equation; problem.append_state = std::bind(&Ephemeris::AppendMasslessBodiesState, _1, std::cref(trajectories)); problem.t_final = t_final; problem.initial_state = &initial_state; AdaptiveStepSize<NewtonianMotionEquation> step_size; step_size.first_time_step = problem.t_final - initial_state.time.value; CHECK_GT(step_size.first_time_step, 0 * Second) << "Flow back to the future: " << problem.t_final << " <= " << initial_state.time.value; step_size.safety_factor = 0.9; step_size.tolerance_to_error_ratio = std::bind(&Ephemeris<Frame>::ToleranceToErrorRatio, std::cref(parameters.length_integration_tolerance_), std::cref(parameters.speed_integration_tolerance_), _1, _2); step_size.max_steps = parameters.max_steps_; auto const outcome = parameters.integrator_->Solve(problem, step_size); // TODO(egg): when we have events in trajectories, we should add a singularity // event at the end if the outcome indicates a singularity // (|VanishingStepSize|). We should not have an event on the trajectory if // |ReachedMaximalStepCount|, since that is not a physical property, but // rather a self-imposed constraint. return outcome == integrators::TerminationCondition::Done && t_final == t; }