SPH
ScriptUtils.h
Go to the documentation of this file.
1 #pragma once
2 
7 
8 #include "post/Analysis.h"
9 #include "quantities/Quantity.h"
10 #include "quantities/Storage.h"
11 #include "system/Factory.h"
12 #include "thread/Scheduler.h"
13 
14 #ifdef SPH_USE_CHAISCRIPT
15 
16 #include <chaiscript/chaiscript.hpp>
17 #include <numeric>
18 
20 
21 namespace Chai {
22 
23 struct Vec3 {
24  Float x, y, z;
25 
26  Vec3() = default;
27 
28  Vec3(const Float x, const Float y, const Float z)
29  : x(x)
30  , y(y)
31  , z(z) {}
32 
33  Vec3(const Vec3& other) {
34  *this = other;
35  }
36 
37  Vec3(const Sph::Vector& v)
38  : x(v[X])
39  , y(v[Y])
40  , z(v[Z]) {}
41 
42  operator Sph::Vector() const {
43  return Sph::Vector(x, y, z);
44  }
45 
46  Vec3& operator=(const Vec3& other) {
47  x = other.x;
48  y = other.y;
49  z = other.z;
50  return *this;
51  }
52 
53  Vec3 operator+(const Vec3& other) const {
54  return { x + other.x, y + other.y, z + other.z };
55  }
56  Vec3 operator-(const Vec3& other) const {
57  return { x - other.x, y - other.y, z - other.z };
58  }
59  Vec3 operator*(const Float f) const {
60  return { x * f, y * f, z * f };
61  }
62  Vec3 operator/(const Float f) const {
63  return { x / f, y / f, z / f };
64  }
65  Vec3& operator+=(const Vec3& other) {
66  *this = *this + other;
67  return *this;
68  }
69 };
70 
71 INLINE Vec3 operator*(const Float f, const Vec3& v) {
72  return v * f;
73 }
74 
75 INLINE Float dot(const Vec3& v1, const Vec3& v2) {
76  return v1.x * v2.x + v1.y * v2.y + v1.z * v2.z;
77 }
78 INLINE Vec3 cross(const Vec3& v1, const Vec3& v2) {
79  return Vec3{ v1.y * v2.z - v1.z * v2.y, v1.z * v2.x - v1.x * v2.z, v1.x * v2.y - v1.y * v2.x };
80 }
81 INLINE Float length(const Vec3& v) {
82  return sqrt(dot(v, v));
83 }
84 INLINE Vec3 min(const Vec3& v1, const Vec3& v2) {
85  return { Sph::min(v1.x, v2.x), Sph::min(v1.z, v2.z), Sph::min(v1.z, v2.z) };
86 }
87 INLINE Vec3 max(const Vec3& v1, const Vec3& v2) {
88  return { Sph::max(v1.x, v2.x), Sph::max(v1.z, v2.z), Sph::max(v1.z, v2.z) };
89 }
90 INLINE Vec3 normalized(const Vec3& v) {
91  return v / length(v);
92 }
93 
94 class Box3 {
95  Vec3 lower, upper;
96 
97 public:
98  Box3()
99  : lower(INFTY, INFTY, INFTY)
100  , upper(-INFTY, -INFTY, -INFTY) {}
101 
102  Box3(const Vec3& lower, const Vec3& upper)
103  : lower(lower)
104  , upper(upper) {}
105 
106  Box3(const Box3& other) {
107  *this = other;
108  }
109 
110  Box3& operator=(const Box3& other) {
111  lower = other.lower;
112  upper = other.upper;
113  return *this;
114  }
115 
116  Vec3 size() const {
117  return upper - lower;
118  }
119 
120  void extend(const Vec3& pos) {
121  lower = min(lower, pos);
122  upper = max(upper, pos);
123  }
124 };
125 
127 class Particles {
128 private:
130  Storage* storage = nullptr;
131 
133  bool owns = true;
134 
136  std::vector<Vec3> positions, velocities, accelerations;
137  std::vector<Float> masses, energies, densities, radii;
138 
139 public:
140  Particles() = default;
141 
142  Particles(const Size particleCnt);
143 
144  ~Particles();
145 
146  void bindToStorage(Storage& input);
147 
148  const Storage& store() const;
149 
150  Size getParticleCnt() const;
151 
152  void resize(const Size particleCnt);
153 
154  std::vector<Float>& getMasses();
155 
156  std::vector<Float>& getEnergies();
157 
158  std::vector<Float>& getDensities();
159 
160  std::vector<Float>& getRadii();
161 
162  std::vector<Vec3>& getPositions();
163 
164  std::vector<Vec3>& getVelocities();
165 
166  std::vector<Vec3>& getAccelerations();
167 
168  Box3 getBoundingBox() const;
169 
170  Float getTotalMass() const;
171 
172  Vec3 getCenterOfMass() const;
173 
174  Vec3 getTotalMomentum() const;
175 
176  Vec3 getTotalAngularMomentum() const;
177 
178  Vec3 getAngularFrequency() const;
179 
180  void merge(Particles& other);
181 };
182 
184 void registerBindings(chaiscript::ChaiScript& chai);
185 
186 } // namespace Chai
187 
189 
190 #endif
Various function for interpretation of the results of a simulation.
NAMESPACE_SPH_BEGIN
Definition: BarnesHut.cpp:13
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
Lut< TValue, TScalar > operator/(const Lut< TValue, TScalar > &lut1, const Lut< TValue, TScalar > &lut2)
Definition: Lut.h:189
constexpr INLINE T max(const T &f1, const T &f2)
Definition: MathBasic.h:20
NAMESPACE_SPH_BEGIN constexpr INLINE T min(const T &f1, const T &f2)
Minimum & Maximum value.
Definition: MathBasic.h:15
INLINE T sqrt(const T f)
Return a squared root of a value.
Definition: MathUtils.h:78
constexpr Float INFTY
Definition: MathUtils.h:38
#define INLINE
Macros for conditional compilation based on selected compiler.
Definition: Object.h:31
#define NAMESPACE_SPH_END
Definition: Object.h:12
Holder of quantity values and their temporal derivatives.
Interface for executing tasks (potentially) asynchronously.
Box getBoundingBox(const Storage &storage, const Float radius)
Convenience function to get the bounding box of all particles.
Definition: Storage.cpp:832
Container for storing particle quantities and materials.
INLINE BasicVector< float > cross(const BasicVector< float > &v1, const BasicVector< float > &v2)
Cross product between two vectors.
Definition: Vector.h:559
BasicVector< Float > Vector
Definition: Vector.h:539
INLINE float dot(const BasicVector< float > &v1, const BasicVector< float > &v2)
Make sure the vector is trivially constructible and destructible, needed for fast initialization of a...
Definition: Vector.h:548
@ Y
Definition: Vector.h:23
@ X
Definition: Vector.h:22
@ Z
Definition: Vector.h:24
Container storing all quantities used within the simulations.
Definition: Storage.h:230
Creating code components based on values from settings.
Difference< TValue1, TValue2 > operator-(const TValue1 &v1, const TValue2 &v2)
Definition: Multipole.h:949
MultiplyByScalar< TValue > operator*(const TValue &v, const Float f)
Definition: Multipole.h:912
Sum< TValue1, TValue2 > operator+(const TValue1 &v1, const TValue2 &v2)
Definition: Multipole.h:933
Vector getAngularFrequency(ArrayView< const Float > m, ArrayView< const Vector > r, ArrayView< const Vector > v, const Vector &r0, const Vector &v0, ArrayView< const Size > idxs=nullptr)
Computes the immediate vector of angular frequency of a rigid body.
Definition: Analysis.cpp:514
Vector getCenterOfMass(ArrayView< const Float > m, ArrayView< const Vector > r, ArrayView< const Size > idxs=nullptr)
Computes the center of mass.
Definition: Analysis.cpp:461