SPH
FragmentationReaccumulation.cpp
Go to the documentation of this file.
1 #include "Sph.h"
2 #include <iostream>
3 
4 using namespace Sph;
5 
7 class LogWriter : public ILogWriter {
8 private:
9  std::string runName;
10 
11 public:
12  LogWriter(const RunSettings& settings)
13  : ILogWriter(makeAuto<StdOutLogger>(), 1._f) {
14  runName = settings.get<std::string>(RunSettingsId::RUN_NAME);
15  }
16 
17  virtual void write(const Storage& UNUSED(storage), const Statistics& stats) override {
18  // Timestep number and current run time
19  const int index = stats.get<int>(StatisticsId::INDEX);
20  const Float time = stats.get<Float>(StatisticsId::RUN_TIME);
21 
22  logger->write(runName, " #", index, " time = ", time);
23  }
24 };
25 
26 // First phase - sets up the colliding bodies and runs the SPH simulation of the fragmentation.
27 class Fragmentation : public IRun {
28 public:
29  virtual void setUp(SharedPtr<Storage> storage) override {
30  settings.set(RunSettingsId::RUN_NAME, std::string("Fragmentation"));
31 
32  // Define the parameters of the target
33  const Float targetRadius = 1.e5_f; // R_pb = 100km;
34  BodySettings targetBody;
35  targetBody
36  .set(BodySettingsId::PARTICLE_COUNT, 10000) // N_pb = 10000
37  .set(BodySettingsId::BODY_RADIUS, targetRadius);
38 
39  // Define the parameters of the impactor
40  const Float impactorRadius = 5.e4_f;
41  BodySettings impactorBody;
42  impactorBody
43  .set(BodySettingsId::BODY_RADIUS, impactorRadius) // R_imp = 50km
44  .set(BodySettingsId::PARTICLE_COUNT, 1250); // N_imp = 1250
45 
46  // Compute position of the impactor from the impact angle
47  const Float impactAngle = 15._f * DEG_TO_RAD;
48  Vector impactorPosition =
49  (targetRadius + impactorRadius) * Vector(cos(impactAngle), sin(impactAngle), 0._f);
50  // Move the impactor to the right, so bodies are not in contact at the start of the simulation
51  impactorPosition[X] += 0.2_f * targetRadius;
52 
53  impactorBody.set(BodySettingsId::BODY_CENTER, impactorPosition);
54 
55  InitialConditions ic(settings);
56  ic.addMonolithicBody(*storage, targetBody);
57 
58  // Using BodyView returned from addMonolithicBody, we can add velocity to the impactor
59  BodyView impactor = ic.addMonolithicBody(*storage, impactorBody);
60  impactor.addVelocity(Vector(-1.e3_f, 0._f, 0._f));
61 
62  // Limit time step by CFL criterion
64 
65  // Run for 100s.
66  settings.set(RunSettingsId::RUN_END_TIME, 100._f);
67 
68  // Save the initial (pre-impact) configuration
69  BinaryOutput io(Path("start.ssf"));
70  Statistics stats;
71  stats.set(StatisticsId::RUN_TIME, 0._f);
72  io.dump(*storage, stats);
73 
74  // Setup custom logging
75  logWriter = makeAuto<LogWriter>(settings);
76  }
77 
78  virtual void tearDown(const Storage& storage, const Statistics& stats) override {
79  // Save the result of the fragmentation phase
80  BinaryOutput io(Path("fragmented.ssf"));
81  io.dump(storage, stats);
82  }
83 };
84 
85 // Second phase - uses the results from the SPH simulation as an input for N-body simulation.
86 class Reaccumulation : public IRun {
87 public:
88  // Here we do not create new particles, the storage already contains particles computed by the
89  // fragmentation phase. Instead we convert the SPH particles to solid spheres (hand-off).
90  virtual void setUp(SharedPtr<Storage> storage) override {
91  settings.set(RunSettingsId::RUN_NAME, std::string("Reaccumulation"));
92 
93  // Create an 'empty' material with default parameters. Note that we do not have to do that in other
94  // examples, as it is done inside InitialConditions object.
95  AutoPtr<IMaterial> material = makeAuto<NullMaterial>(BodySettings::getDefaults());
96 
97  // Create a new (empty) particle storage with our material.
98  Storage spheres(std::move(material));
99 
100  // Get positions and smoothing lengths of SPH particles;
101  // the smoothing lengths can be accessed using r_sph[i][H].
103 
104  // Get velocities of SPH particles (= derivatives of positions)
105  Array<Vector>& v_sph = storage->getDt<Vector>(QuantityId::POSITION);
106 
107  // Get masses and densities of SPH particles
108  Array<Float>& m_sph = storage->getValue<Float>(QuantityId::MASS);
109  Array<Float>& rho_sph = storage->getValue<Float>(QuantityId::DENSITY);
110 
111 
112  // Insert the positions of spheres into the new particle storage. We also want to initialize
113  // velocities and accelerations of the sphere, hence OrderEnum::SECOND.
114  // Note that Array is non-copyable, so we have to clone the data.
116 
117  // Copy also the velocities. Since the velocities were already initialized above (to zeros),
118  // we can simply override the array.
119  spheres.getDt<Vector>(QuantityId::POSITION) = v_sph.clone();
120 
121  // Insert the masses. Masses of spheres correspond to masses of SPH particles, so just copy them.
122  // Note that masses are constant during the simulation (there are no mass derivatives), hence
123  // OrderEnum::ZERO.
124  spheres.insert<Float>(QuantityId::MASS, OrderEnum::ZERO, m_sph.clone());
125 
126  // In N-body simulation, we also need the angular frequencies of particles - initialize them to zero.
128 
129  // Finally, we need to set the radii of the created spheres. Let's choose the radii so that the
130  // volume of the sphere is the same as the volume of the SPH particles.
131  Array<Vector>& r_nbody = spheres.getValue<Vector>(QuantityId::POSITION);
132  for (Size i = 0; i < r_nbody.size(); ++i) {
133  r_nbody[i][H] = Sph::cbrt(3._f * m_sph[i] / (4._f * PI * rho_sph[i]));
134  }
135 
136  // Once all quantities needed for N-body simulation are set up, replace the particles originally
137  // passed into the function; from this point, we do not need SPH data anymore, so we can deallocate
138  // them.
139  *storage = std::move(spheres);
140 
141 
142  // Here we need to select a non-default solver - instead of SPH solver, use N-body solver
146 
147  solver = makeAuto<HardSphereSolver>(*scheduler, settings);
148 
149  // For manually created solvers, it is necessary to call function create for every material in the
150  // storage. In this case, we only have one "sort" of particles, so we call the function just once for
151  // material with index 0.
152  solver->create(*storage, storage->getMaterial(0));
153 
154  // Use the leaf-frog integrator
156 
157  // Limit the time step by accelerations. By default, the time step is also limited by CFL criterion,
158  // which requires sound speed computed for every particle. Here, we do not store the sound speed, so
159  // the code would report a runtime error if we used the CFL criterion.
162 
163  // Set output quantities
167 
168  // Run for 6 hours
169  settings.set(RunSettingsId::RUN_END_TIME, 60._f * 60._f * 6._f);
170 
171  // Setup custom logging
172  logWriter = makeAuto<LogWriter>(settings);
173  }
174 
175  virtual void tearDown(const Storage& storage, const Statistics& stats) override {
176  // Save the result of the reaccumulation phase
177  BinaryOutput io(Path("reaccumulated.ssf"));
178  io.dump(storage, stats);
179  }
180 };
181 
182 int main() {
183  try {
184  Storage storage;
185 
186  // Set up and run the fragmentation phase
187  Fragmentation fragmentation;
188  fragmentation.run(storage);
189 
190  // Now storage contains SPH particles, we can pass them as an input to to the reaccumulation phase
191  Reaccumulation reaccumulation;
192  reaccumulation.run(storage);
193 
194  } catch (const Exception& e) {
195  std::cout << "Error during simulation: " << e.what() << std::endl;
196  return -1;
197  }
198  return 0;
199 }
INLINE AutoPtr< T > makeAuto(TArgs &&... args)
Definition: AutoPtr.h:124
uint32_t Size
Integral type used to index arrays (by default).
Definition: Globals.h:16
double Float
Precision used withing the code. Use Float instead of float or double where precision is important.
Definition: Globals.h:13
INLINE Float cbrt(const Float f)
Returns a cubed root of a value.
Definition: MathUtils.h:84
INLINE T sin(const T f)
Definition: MathUtils.h:296
INLINE T cos(const T f)
Definition: MathUtils.h:291
constexpr Float DEG_TO_RAD
Definition: MathUtils.h:369
constexpr Float PI
Mathematical constants.
Definition: MathUtils.h:361
#define UNUSED(x)
Definition: Object.h:37
@ VELOCITY
Current velocities of particles, always a vector quantity.
@ POSITION
Positions of particles, always a vector quantity.
@ SMOOTHING_LENGTH
Smoothing lenghts of particles.
@ INDEX
Index of particle.
@ MASS
Particle masses, always a scalar quantity.
@ POSITION
Positions (velocities, accelerations) of particles, always a vector quantity,.
@ DENSITY
Density, always a scalar quantity.
@ MASS
Paricles masses, always a scalar quantity.
@ SECOND
Quantity with 1st and 2nd derivative.
@ ZERO
Quantity without derivatives, or "zero order" of quantity.
Includes common headers.
@ RUN_TIME
Current time of the simulation in code units. Does not necessarily have to be 0 when run starts.
@ INDEX
Current number of time step, indexed from 0.
BasicVector< Float > Vector
Definition: Vector.h:539
@ H
Definition: Vector.h:25
@ X
Definition: Vector.h:22
Generic dynamically allocated resizable storage.
Definition: Array.h:43
Array clone() const
Performs a deep copy of all elements of the array.
Definition: Array.h:147
Output saving data to binary data without loss of precision.
Definition: Output.h:335
virtual Expected< Path > dump(const Storage &storage, const Statistics &stats) override
Saves data from particle storage into the file.
Definition: Output.cpp:512
Non-owning view of particles belonging to the same body.
Definition: Initial.h:20
BodyView & addVelocity(const Vector &v)
Adds a velocity vector to all particles of the body.
Definition: Initial.cpp:38
Generic exception.
Definition: Exceptions.h:10
virtual const char * what() const noexcept
Definition: Exceptions.h:18
virtual void setUp(SharedPtr< Storage > storage) override
Prepares the run, creates logger, output, ...
virtual void tearDown(const Storage &storage, const Statistics &stats) override
Called after the run.
Base class for objects logging run statistics.
Definition: LogWriter.h:10
Defines the interface for a run.
Definition: IRun.h:61
Statistics run(Storage &storage)
Runs the simulation.
Definition: IRun.cpp:187
Object for adding one or more bodies with given material into a Storage.
Definition: Initial.h:106
BodyView addMonolithicBody(Storage &storage, const BodySettings &body)
Creates a monolithic body by filling given domain with particles.
Definition: Initial.cpp:81
Custom writer of simulation log.
virtual void write(const Storage &UNUSED(storage), const Statistics &stats) override
LogWriter(const RunSettings &settings)
Object representing a path on a filesystem.
Definition: Path.h:17
virtual void setUp(SharedPtr< Storage > storage) override
Prepares the run, creates logger, output, ...
virtual void tearDown(const Storage &storage, const Statistics &stats) override
Called after the run.
TValue get(const TEnum idx, std::enable_if_t<!std::is_enum< std::decay_t< TValue >>::value, int >=0) const
Returns a value of given type from the settings.
Definition: Settings.h:326
static const Settings & getDefaults()
\brief Returns a reference to object containing default values of all settings.
Settings & set(const TEnum idx, TValue &&value, std::enable_if_t<!std::is_enum< std::decay_t< TValue >>::value, int >=0)
Saves a value into the settings.
Definition: Settings.h:226
Object holding various statistics about current run.
Definition: Statistics.h:22
TValue get(const StatisticsId idx) const
Returns value of a statistic.
Definition: Statistics.h:88
Statistics & set(const StatisticsId idx, TValue &&value)
Sets new values of a statistic.
Definition: Statistics.h:52
Standard output logger.
Definition: Logger.h:138
Container storing all quantities used within the simulations.
Definition: Storage.h:230
Quantity & insert(const QuantityId key, const OrderEnum order, const TValue &defaultValue)
Creates a quantity in the storage, given its key, value type and order.
Definition: Storage.cpp:270
Array< TValue > & getDt(const QuantityId key)
Retrieves a quantity derivative from the storage, given its key and value type.
Definition: Storage.cpp:217
MaterialView getMaterial(const Size matIdx) const
Returns an object containing a reference to given material.
Definition: Storage.cpp:366
Array< TValue > & getValue(const QuantityId key)
Retrieves a quantity values from the storage, given its key and value type.
Definition: Storage.cpp:191
@ COURANT
Time step determined using CFL condition.
@ ACCELERATION
Time step computed from ratio of acceleration and smoothing length.
@ LEAP_FROG
Leap-frog 2nd-order integration.
@ BODY_CENTER
Center point of the body. Currently used only by StabilizationSolver.
@ PARTICLE_COUNT
Number of SPH particles in the body.
@ RUN_OUTPUT_QUANTITIES
List of quantities to write to text output. Binary output always stores all quantitites.
@ TIMESTEPPING_MAX_TIMESTEP
@ COLLISION_HANDLER
Specifies how the collisions of particles should be handler; see CollisionHandlerEnum.
@ TIMESTEPPING_INTEGRATOR
Selected timestepping integrator.
@ COLLISION_OVERLAP
Specifies how particle overlaps should be handled.
@ RUN_NAME
User-specified name of the run, used in some output files.
@ GRAVITY_KERNEL
Gravity smoothing kernel.
Definition: MemoryPool.h:5