MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::navier::NavierParticles Class Reference

Transient Navier-Stokes fluid-particles solver. More...

#include <navier_particles.hpp>

Collaboration diagram for mfem::navier::NavierParticles:
[legend]

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.
 
ParticleSetGetParticles ()
 Get reference to the active ParticleSet.
 
ParticleSetGetInactiveParticles ()
 Get reference to the inactive ParticleSet.
 
ParticleVectorKappa ()
 Get reference to the κ ParticleVector.
 
ParticleVectorZeta ()
 Get reference to the ζ ParticleVector.
 
ParticleVectorGamma ()
 Get reference to the γ ParticleVector.
 
ParticleVectorU (int nm=0)
 Get reference to the fluid velocity-interpolated ParticleVector at time n - nm .
 
ParticleVectorV (int nm=0)
 Get reference to the particle velocity ParticleVector at time n - nm .
 
ParticleVectorW (int nm=0)
 Get reference to the fluid vorticity-interpolated ParticleVector at time n - nm .
 
ParticleVectorX (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_tdthist {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< BCVariantbcs
 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)
 

Detailed Description

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.

Member Typedef Documentation

◆ BCVariant

Definition at line 113 of file navier_particles.hpp.

Constructor & Destructor Documentation

◆ NavierParticles()

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.

Member Function Documentation

◆ Add2DRecirculationBC()

void mfem::navier::NavierParticles::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 )
inline

Add a 2D recirculation / one-way periodic boundary condition.

Warning
Both normals must be facing into the domain. See Get2DNormal() for details on the normal direction.
Parameters
[in]inlet_startInlet line segment start point.
[in]inlet_endInlet line segment end point.
[in]invert_inlet_normalInvert direction of the inlet normal.
[in]outlet_startOutlet line segment start point.
[in]outlet_endOutlet line segment end point.
[in]invert_outlet_normalInvert direction of the outlet normal.

Definition at line 297 of file navier_particles.hpp.

◆ Add2DReflectionBC()

void mfem::navier::NavierParticles::Add2DReflectionBC ( const Vector & line_start,
const Vector & line_end,
real_t e,
bool invert_normal )
inline

Add a 2D wall reflective boundary condition.

Parameters
[in]line_startWall line segment start point.
[in]line_endWall line segment end point.
[in]eBoundary collision reconstitution constant. 1 for perfectly elastic.
[in]invert_normalTrue if left normal points out of domain. False if left normal points into domain.

Definition at line 280 of file navier_particles.hpp.

◆ Apply2DRecirculationBC()

void mfem::navier::NavierParticles::Apply2DRecirculationBC ( const RecirculationBC_2D & bc)
protected

Apply 2D recirculation BCs.

Definition at line 243 of file navier_particles.cpp.

◆ Apply2DReflectionBC()

void mfem::navier::NavierParticles::Apply2DReflectionBC ( const ReflectionBC_2D & bc)
protected

Apply 2D reflection BCs to update particle positions and velocities.

Definition at line 200 of file navier_particles.cpp.

◆ ApplyBCs()

void mfem::navier::NavierParticles::ApplyBCs ( )
protected

Apply all BCs in bcs.

Definition at line 298 of file navier_particles.cpp.

◆ DeactivateLostParticles()

void mfem::navier::NavierParticles::DeactivateLostParticles ( bool findpts)
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.

Parameters
[in]findptsIf true, call FindPointsGSLIB::FindPoints prior to deactivation (if particle coordinates out of sync with FindPointsGSLIB)

Definition at line 435 of file navier_particles.cpp.

◆ Gamma()

ParticleVector & mfem::navier::NavierParticles::Gamma ( )
inline

Get reference to the γ ParticleVector.

Definition at line 229 of file navier_particles.hpp.

◆ Get2DNormal()

void mfem::navier::NavierParticles::Get2DNormal ( const Vector & p1,
const Vector & p2,
bool inv_normal,
Vector & normal )
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.

◆ Get2DSegmentIntersection()

bool mfem::navier::NavierParticles::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 )
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 .

Parameters
[in]s1_startFirst line segment start point.
[in]s1_endFirst line segment end point.
[in]s2_startSecond line segment start point.
[in]s2_endSecond line segment end point.
[in]x_intComputed intersection point (if it exists)
[in]t1_ptr(Optional) Computed t1
[in]t2_ptr(Optional) Computed t2
Returns
true if intersection point exists, false otherwise

Definition at line 155 of file navier_particles.cpp.

◆ GetInactiveParticles()

ParticleSet & mfem::navier::NavierParticles::GetInactiveParticles ( )
inline

Get reference to the inactive ParticleSet.

Definition at line 214 of file navier_particles.hpp.

◆ GetParticles()

ParticleSet & mfem::navier::NavierParticles::GetParticles ( )
inline

Get reference to the active ParticleSet.

Definition at line 211 of file navier_particles.hpp.

◆ InterpolateUW()

void mfem::navier::NavierParticles::InterpolateUW ( const ParGridFunction & u_gf,
const ParGridFunction & w_gf )

Interpolate fluid velocity and vorticity onto current particles' location.

Parameters
[in]u_gfFluid velocity on fluid mesh.
[in]w_gfFluid vorticity on fluid mesh.

Definition at line 421 of file navier_particles.cpp.

◆ Kappa()

ParticleVector & mfem::navier::NavierParticles::Kappa ( )
inline

Get reference to the κ ParticleVector.

Definition at line 217 of file navier_particles.hpp.

◆ Order()

Array< int > & mfem::navier::NavierParticles::Order ( )
inline

Get reference to the order Array<int>.

Definition at line 268 of file navier_particles.hpp.

◆ ParticleStep2D()

void mfem::navier::NavierParticles::ParticleStep2D ( const real_t & dt,
int p )
protected

2D particle step for particle at index p

Definition at line 74 of file navier_particles.cpp.

◆ SetTimeIntegrationCoefficients()

void mfem::navier::NavierParticles::SetTimeIntegrationCoefficients ( )
protected

Update beta_k and alpha_k using current dthist.

Definition at line 23 of file navier_particles.cpp.

◆ Setup()

void mfem::navier::NavierParticles::Setup ( const real_t dt)
inline

Set initial timestep in time history array.

Definition at line 191 of file navier_particles.hpp.

◆ Step()

void mfem::navier::NavierParticles::Step ( const real_t dt,
const ParGridFunction & u_gf,
const ParGridFunction & w_gf )

Step the particles in time.

Parameters
[in]dt
[in]u_gfFluid velocity on fluid mesh.
[in]w_gfFluid vorticity on fluid mesh.

Definition at line 369 of file navier_particles.cpp.

◆ U()

ParticleVector & mfem::navier::NavierParticles::U ( int nm = 0)
inline

Get reference to the fluid velocity-interpolated ParticleVector at time n - nm .

Definition at line 237 of file navier_particles.hpp.

◆ V()

ParticleVector & mfem::navier::NavierParticles::V ( int nm = 0)
inline

Get reference to the particle velocity ParticleVector at time n - nm .

Definition at line 244 of file navier_particles.hpp.

◆ W()

ParticleVector & mfem::navier::NavierParticles::W ( int nm = 0)
inline

Get reference to the fluid vorticity-interpolated ParticleVector at time n - nm .

Definition at line 253 of file navier_particles.hpp.

◆ X()

ParticleVector & mfem::navier::NavierParticles::X ( int nm = 0)
inline

Get reference to the position ParticleVector at time n - nm .

Definition at line 260 of file navier_particles.hpp.

◆ Zeta()

ParticleVector & mfem::navier::NavierParticles::Zeta ( )
inline

Get reference to the ζ ParticleVector.

Definition at line 223 of file navier_particles.hpp.

Member Data Documentation

◆ alpha_k

std::array<Array<real_t>, 3> mfem::navier::NavierParticles::alpha_k
protected

EXTk coefficients, k=1,2,3.

Definition at line 69 of file navier_particles.hpp.

◆ bcs

std::vector<BCVariant> mfem::navier::NavierParticles::bcs
protected

std::vector of all boundary condition structs

Definition at line 116 of file navier_particles.hpp.

◆ beta_k

std::array<Array<real_t>, 3> mfem::navier::NavierParticles::beta_k
protected

BDFk coefficients, k=1,2,3.

Definition at line 66 of file navier_particles.hpp.

◆ C

Vector mfem::navier::NavierParticles::C
protected

Definition at line 183 of file navier_particles.hpp.

◆ dthist

Array<real_t> mfem::navier::NavierParticles::dthist {0.0, 0.0, 0.0}
protected

Timestep history.

Definition at line 63 of file navier_particles.hpp.

◆ finder

FindPointsGSLIB mfem::navier::NavierParticles::finder
protected

Definition at line 60 of file navier_particles.hpp.

◆ fluid_particles

ParticleSet mfem::navier::NavierParticles::fluid_particles
protected

Active fluid particle set.

Definition at line 54 of file navier_particles.hpp.

◆ fp_idx

struct mfem::navier::NavierParticles::FluidParticleIndices mfem::navier::NavierParticles::fp_idx
protected

◆ inactive_fluid_particles

ParticleSet mfem::navier::NavierParticles::inactive_fluid_particles
protected

Inactive fluid particle set. Particles that leave the domain are added here.

Definition at line 58 of file navier_particles.hpp.

◆ N_HIST

int mfem::navier::NavierParticles::N_HIST = 4
staticconstexprprotected

Number of timesteps to store (includes current)

Definition at line 72 of file navier_particles.hpp.

◆ r

Vector mfem::navier::NavierParticles::r
mutableprotected

Definition at line 183 of file navier_particles.hpp.

◆ up

Vector mfem::navier::NavierParticles::up
mutableprotected

Definition at line 182 of file navier_particles.hpp.

◆ vp

Vector mfem::navier::NavierParticles::vp
protected

Definition at line 182 of file navier_particles.hpp.

◆ xp

Vector mfem::navier::NavierParticles::xp
protected

Definition at line 182 of file navier_particles.hpp.

◆ xpn

Vector mfem::navier::NavierParticles::xpn
protected

Definition at line 182 of file navier_particles.hpp.


The documentation for this class was generated from the following files: