![]() |
MFEM v4.9.0
Finite element discretization library
|
Transient Navier-Stokes fluid-particles solver. More...
#include <navier_particles.hpp>
Classes | |
| struct | FluidParticleIndices |
| struct | RecirculationBC_2D |
| 2D recirculation boundary condition struct More... | |
| struct | ReflectionBC_2D |
| 2D wall reflection boundary condition struct More... | |
Public Member Functions | |
| NavierParticles (MPI_Comm comm, int num_particles, Mesh &m) | |
Initialize NavierParticles with num_particles using fluid mesh m . | |
| void | Setup (const real_t dt) |
| Set initial timestep in time history array. | |
| void | Step (const real_t dt, const ParGridFunction &u_gf, const ParGridFunction &w_gf) |
| Step the particles in time. | |
| void | InterpolateUW (const ParGridFunction &u_gf, const ParGridFunction &w_gf) |
| Interpolate fluid velocity and vorticity onto current particles' location. | |
| ParticleSet & | GetParticles () |
| Get reference to the active ParticleSet. | |
| ParticleSet & | GetInactiveParticles () |
| Get reference to the inactive ParticleSet. | |
| ParticleVector & | Kappa () |
| Get reference to the κ ParticleVector. | |
| ParticleVector & | Zeta () |
| Get reference to the ζ ParticleVector. | |
| ParticleVector & | Gamma () |
| Get reference to the γ ParticleVector. | |
| ParticleVector & | U (int nm=0) |
Get reference to the fluid velocity-interpolated ParticleVector at time n - nm . | |
| ParticleVector & | V (int nm=0) |
Get reference to the particle velocity ParticleVector at time n - nm . | |
| ParticleVector & | W (int nm=0) |
Get reference to the fluid vorticity-interpolated ParticleVector at time n - nm . | |
| ParticleVector & | X (int nm=0) |
Get reference to the position ParticleVector at time n - nm . | |
| Array< int > & | Order () |
| Get reference to the order Array<int>. | |
| void | Add2DReflectionBC (const Vector &line_start, const Vector &line_end, real_t e, bool invert_normal) |
| Add a 2D wall reflective boundary condition. | |
| void | Add2DRecirculationBC (const Vector &inlet_start, const Vector &inlet_end, bool invert_inlet_normal, const Vector &outlet_start, const Vector &outlet_end, bool invert_outlet_normal) |
| Add a 2D recirculation / one-way periodic boundary condition. | |
Protected Types | |
| using | BCVariant = std::variant<ReflectionBC_2D, RecirculationBC_2D> |
Protected Member Functions | |
| void | SetTimeIntegrationCoefficients () |
| Update beta_k and alpha_k using current dthist. | |
| void | ParticleStep2D (const real_t &dt, int p) |
2D particle step for particle at index p | |
| void | Apply2DReflectionBC (const ReflectionBC_2D &bc) |
| Apply 2D reflection BCs to update particle positions and velocities. | |
| void | Apply2DRecirculationBC (const RecirculationBC_2D &bc) |
| Apply 2D recirculation BCs. | |
| void | ApplyBCs () |
| Apply all BCs in bcs. | |
| void | DeactivateLostParticles (bool findpts) |
| Move any particles that have left the domain to the inactive ParticleSet. | |
Static Protected Member Functions | |
| static void | Get2DNormal (const Vector &p1, const Vector &p2, bool inv_normal, Vector &normal) |
| Given two 2D points, get the unit normal to the line connecting them. | |
| static bool | Get2DSegmentIntersection (const Vector &s1_start, const Vector &s1_end, const Vector &s2_start, const Vector &s2_end, Vector &x_int, real_t *t1_ptr=nullptr, real_t *t2_ptr=nullptr) |
| Given two 2D segments, get the intersection point if it exists. | |
Protected Attributes | |
| ParticleSet | fluid_particles |
| Active fluid particle set. | |
| ParticleSet | inactive_fluid_particles |
| FindPointsGSLIB | finder |
| Array< real_t > | dthist {0.0, 0.0, 0.0} |
| Timestep history. | |
| std::array< Array< real_t >, 3 > | beta_k |
| BDFk coefficients, k=1,2,3. | |
| std::array< Array< real_t >, 3 > | alpha_k |
| EXTk coefficients, k=1,2,3. | |
| struct mfem::navier::NavierParticles::FluidParticleIndices | fp_idx |
| std::vector< BCVariant > | bcs |
| std::vector of all boundary condition structs | |
| Vector | up |
| Vector | vp |
| Vector | xpn |
| Vector | xp |
| Vector | r |
| Vector | C |
Static Protected Attributes | |
| static constexpr int | N_HIST = 4 |
| Number of timesteps to store (includes current) | |
Transient Navier-Stokes fluid-particles solver.
This class implements a Lagrangian point particle tracking model from Dutta, Som (2017). Ph.D. Dissertation: Bulle-Effect and its implications for morphodynamics of river diversions. https://www.ideals.illinois.edu/items/102343.
In short, the particles are advanced in time by solving the ODE: dv/dt = κ(u-v) - γ ê + ζ (u-v) x ω, dx/dt = v, where x and v are the particle location and velocity, u is the fluid velocity at the particle location, ω is the fluid vorticity at the particle location, κ depends on the drag properties of the particle, ζ depends on the lift properties, and γ and ê depend on body forces such as gravity.
The model from Dutta et al. is general but this implementation is currently limited to 2D problems. Simple reflection and recirculation boundary conditions are also supported.
Definition at line 50 of file navier_particles.hpp.
|
protected |
Definition at line 113 of file navier_particles.hpp.
| mfem::navier::NavierParticles::NavierParticles | ( | MPI_Comm | comm, |
| int | num_particles, | ||
| Mesh & | m ) |
Initialize NavierParticles with num_particles using fluid mesh m .
Definition at line 319 of file navier_particles.cpp.
|
inline |
Add a 2D recirculation / one-way periodic boundary condition.
| [in] | inlet_start | Inlet line segment start point. |
| [in] | inlet_end | Inlet line segment end point. |
| [in] | invert_inlet_normal | Invert direction of the inlet normal. |
| [in] | outlet_start | Outlet line segment start point. |
| [in] | outlet_end | Outlet line segment end point. |
| [in] | invert_outlet_normal | Invert direction of the outlet normal. |
Definition at line 297 of file navier_particles.hpp.
|
inline |
Add a 2D wall reflective boundary condition.
| [in] | line_start | Wall line segment start point. |
| [in] | line_end | Wall line segment end point. |
| [in] | e | Boundary collision reconstitution constant. 1 for perfectly elastic. |
| [in] | invert_normal | True if left normal points out of domain. False if left normal points into domain. |
Definition at line 280 of file navier_particles.hpp.
|
protected |
Apply 2D recirculation BCs.
Definition at line 243 of file navier_particles.cpp.
|
protected |
Apply 2D reflection BCs to update particle positions and velocities.
Definition at line 200 of file navier_particles.cpp.
|
protected |
Apply all BCs in bcs.
Definition at line 298 of file navier_particles.cpp.
|
protected |
Move any particles that have left the domain to the inactive ParticleSet.
This method uses the FindPointsGSLIB object internal to the class to detect if particles are within the domain or not.
| [in] | findpts | If true, call FindPointsGSLIB::FindPoints prior to deactivation (if particle coordinates out of sync with FindPointsGSLIB) |
Definition at line 435 of file navier_particles.cpp.
|
inline |
Get reference to the γ ParticleVector.
Definition at line 229 of file navier_particles.hpp.
|
staticprotected |
Given two 2D points, get the unit normal to the line connecting them.
For inv_normal false, the normal is 90 degrees (CCW) from the line connecting p1 to p2 . For inv_normal true, the normal is 90 degrees (CW) from the line connecting p1 to p2 .
Definition at line 136 of file navier_particles.cpp.
|
staticprotected |
Given two 2D segments, get the intersection point if it exists.
Specifically this function solves the following system of equations: r_1 = s1_start + t1*[s1_end - s1_start] r_2 = s2_start + t2*[s2_end - s2_start]
and then checks if 0<=t1<=1 and 0<=t2<=1 .
| [in] | s1_start | First line segment start point. |
| [in] | s1_end | First line segment end point. |
| [in] | s2_start | Second line segment start point. |
| [in] | s2_end | Second line segment end point. |
| [in] | x_int | Computed intersection point (if it exists) |
| [in] | t1_ptr | (Optional) Computed t1 |
| [in] | t2_ptr | (Optional) Computed t2 |
Definition at line 155 of file navier_particles.cpp.
|
inline |
Get reference to the inactive ParticleSet.
Definition at line 214 of file navier_particles.hpp.
|
inline |
Get reference to the active ParticleSet.
Definition at line 211 of file navier_particles.hpp.
| void mfem::navier::NavierParticles::InterpolateUW | ( | const ParGridFunction & | u_gf, |
| const ParGridFunction & | w_gf ) |
Interpolate fluid velocity and vorticity onto current particles' location.
| [in] | u_gf | Fluid velocity on fluid mesh. |
| [in] | w_gf | Fluid vorticity on fluid mesh. |
Definition at line 421 of file navier_particles.cpp.
|
inline |
Get reference to the κ ParticleVector.
Definition at line 217 of file navier_particles.hpp.
|
inline |
Get reference to the order Array<int>.
Definition at line 268 of file navier_particles.hpp.
|
protected |
2D particle step for particle at index p
Definition at line 74 of file navier_particles.cpp.
|
protected |
Update beta_k and alpha_k using current dthist.
Definition at line 23 of file navier_particles.cpp.
|
inline |
Set initial timestep in time history array.
Definition at line 191 of file navier_particles.hpp.
| void mfem::navier::NavierParticles::Step | ( | const real_t | dt, |
| const ParGridFunction & | u_gf, | ||
| const ParGridFunction & | w_gf ) |
Step the particles in time.
| [in] | dt | |
| [in] | u_gf | Fluid velocity on fluid mesh. |
| [in] | w_gf | Fluid vorticity on fluid mesh. |
Definition at line 369 of file navier_particles.cpp.
|
inline |
Get reference to the fluid velocity-interpolated ParticleVector at time n - nm .
Definition at line 237 of file navier_particles.hpp.
|
inline |
Get reference to the particle velocity ParticleVector at time n - nm .
Definition at line 244 of file navier_particles.hpp.
|
inline |
Get reference to the fluid vorticity-interpolated ParticleVector at time n - nm .
Definition at line 253 of file navier_particles.hpp.
|
inline |
Get reference to the position ParticleVector at time n - nm .
Definition at line 260 of file navier_particles.hpp.
|
inline |
Get reference to the ζ ParticleVector.
Definition at line 223 of file navier_particles.hpp.
EXTk coefficients, k=1,2,3.
Definition at line 69 of file navier_particles.hpp.
|
protected |
std::vector of all boundary condition structs
Definition at line 116 of file navier_particles.hpp.
BDFk coefficients, k=1,2,3.
Definition at line 66 of file navier_particles.hpp.
|
protected |
Definition at line 183 of file navier_particles.hpp.
Timestep history.
Definition at line 63 of file navier_particles.hpp.
|
protected |
Definition at line 60 of file navier_particles.hpp.
|
protected |
Active fluid particle set.
Definition at line 54 of file navier_particles.hpp.
|
protected |
|
protected |
Inactive fluid particle set. Particles that leave the domain are added here.
Definition at line 58 of file navier_particles.hpp.
|
staticconstexprprotected |
Number of timesteps to store (includes current)
Definition at line 72 of file navier_particles.hpp.
|
mutableprotected |
Definition at line 183 of file navier_particles.hpp.
|
mutableprotected |
Definition at line 182 of file navier_particles.hpp.
|
protected |
Definition at line 182 of file navier_particles.hpp.
|
protected |
Definition at line 182 of file navier_particles.hpp.
|
protected |
Definition at line 182 of file navier_particles.hpp.