Template Class Inclusions¶
Defined in File inclusions.h
Inheritance Relationships¶
Base Type¶
public ParameterAcceptor
Class Documentation¶
-
template<int spacedim>
class Inclusions : public ParameterAcceptor¶ Class for handling inclusions in an immersed boundary method.
This class provides functionality for handling inclusions in an immersed boundary method. It stores the positions, radii, and Fourier coefficients of the inclusions, and provides methods for computing their normals, center, and direction. It also provides methods for initializing the inclusions from a file, setting up particles for the inclusions, and getting degrees of freedom associated with each inclusion.
Public Functions
-
inline Inclusions(const unsigned int n_vector_components = 1)¶
Class for computing the inclusions of a given mesh.
- Parameters:
n_vector_components – Number of vector components.
-
inline types::global_dof_index n_dofs() const¶
Returns the number of degrees of freedom in the system.
- Returns:
The number of degrees of freedom.
-
inline types::global_dof_index n_particles() const¶
Returns the number of particles in the system.
- Returns:
The number of particles in the system.
-
inline types::global_dof_index n_inclusions() const¶
Returns the number of inclusions in the mesh.
- Returns:
The number of inclusions in the mesh.
-
inline unsigned int n_global_segments() const¶
-
inline unsigned int n_local_segments() const¶
-
inline unsigned int n_dofs_per_inclusion() const¶
Number of degrees of freedom associated to each inclusion.
- Returns:
unsigned int
-
inline void initialize()¶
Reinit all reference points and normals, and read inclusions file.
-
inline std::vector<types::global_dof_index> get_dof_indices(const types::global_dof_index &quadrature_id) const¶
Returns the degrees of freedom indices associated with a given quadrature point.
- Parameters:
quadrature_id – The global index of the quadrature point.
- Returns:
A vector containing the degrees of freedom indices.
-
inline void setup_inclusions_particles(const parallel::distributed::Triangulation<spacedim> &tria)¶
Sets up the inclusions particles for the given triangulation.
- Parameters:
tria – The triangulation to set up the inclusions particles for.
-
inline types::global_dof_index get_inclusion_id(const types::global_dof_index &quadrature_id) const¶
Quadrature id to inclusion id.
- Parameters:
quadrature_id –
- Returns:
const types::global_dof_index
-
inline unsigned int get_component(const types::global_dof_index &dof_index) const¶
Get the ith component for the given dof index.
- Parameters:
dof_index – A number in [0,n_dofs())
- Returns:
unsigned int The index of the current component
-
inline unsigned int get_fourier_component(const types::global_dof_index &dof_index) const¶
Get the ith Fourier component for the given dof index.
- Parameters:
dof_index – A number in [0,n_dofs())
- Returns:
unsigned int The index of the current component
-
inline double get_inclusion_data(const types::global_dof_index &inclusion_id, const types::global_dof_index &dof_index) const¶
Get the Fourier data for the given local dof index.
- Parameters:
inclusion_id – A number in [0,n_inclusions())
dof_index – A number in [0,n_dofs_per_inclusion())
- Returns:
unsigned int The index of the current component
-
inline double get_inclusion_data(const types::global_dof_index &inclusion_id, const types::global_dof_index &dof_index, const Point<spacedim> &point) const¶
Get the Fourier data for the given local dof index, computed at the given quadrature point.
- Parameters:
inclusion_id – A number in [0,n_inclusions())
dof_index – A number in [0,n_dofs_per_inclusion())
point – The inclusion quadrature point location
- Returns:
unsigned int The index of the current component
-
inline const Tensor<1, spacedim> &get_normal(const types::global_dof_index &quadrature_id) const¶
Get the normal.
- Parameters:
quadrature_id –
- Returns:
const Tensor<1, spacedim>&
-
inline double get_JxW(const types::global_dof_index &particle_id) const¶
return weight for integration normal direction of an inclusion at a quadrature point
- Parameters:
particle_id –
- Returns:
const Tensor<1, spacedim>&
-
inline const std::vector<double> &get_fe_values(const types::global_dof_index particle_id) const¶
Get a list of fe values.
- Parameters:
particle_id –
- Returns:
const std::vector<double>&
-
inline Point<spacedim> get_center(const types::global_dof_index &inclusion_id) const¶
Get the center of the inclusion.
- Parameters:
inclusion_id –
- Returns:
Point<spacedim>
-
inline double get_radius(const types::global_dof_index &inclusion_id) const¶
Get the radius of the inclusion.
- Parameters:
inclusion_id –
- Returns:
double
-
inline double get_section_measure(const types::global_dof_index &inclusion_id) const¶
Get the measure of the section of the inclusion.
- Parameters:
inclusion_id –
- Returns:
double
-
inline Tensor<1, spacedim> get_direction(const types::global_dof_index &inclusion_id) const¶
Get the direction of the inclusion.
- Parameters:
inclusion_id –
- Returns:
Tensor<1, spacedim>
-
inline Tensor<2, spacedim> get_rotation(const types::global_dof_index &inclusion_id) const¶
Get the rotation of the inclusion.
- Parameters:
inclusion_id –
- Returns:
Tensor<2, spacedim>
-
inline const std::vector<Point<spacedim>> &get_current_support_points(const types::global_dof_index &inclusion_id) const¶
Return support points of one inclusion in current configuration.
-
inline void output_particles(const std::string &filename) const¶
print the inclusions in parallel on a .vtu file
- Parameters:
filename –
-
inline void update_displacement_hdf5()¶
Update the displacement data after the initialization reading from a hdf5 file.
-
inline void update_inclusions_data(std::vector<double> new_data)¶
Update the displacement data after the initialization, 3D case only.
- Parameters:
new_data – vector of lenght equal to the number of inclusions or to the number of vessels, elements of the vector are the new values to be assigned for normal expansion
-
inline void update_inclusions_data(std::vector<std::vector<double>> new_data)¶
Update inclusion data by vessel-wise sampled profiles.
-
inline int get_vesselID(const types::global_dof_index &inclusion_id) const¶
3D return the vessel that a given inclusion belongs to, 2D return 0
- Parameters:
inclusion_id –
- Returns:
unsigned int of vessel id
-
inline void update_single_inclusion_data_along_normal(const types::global_dof_index &inclusion_id, const double nd)¶
Update one inclusion using a scalar normal-displacement value.
-
inline void update_single_vessel_data(const types::global_dof_index&, const std::vector<double>)¶
Placeholder for vessel-wise updates with vector-valued payloads.
-
inline unsigned int get_n_vessels() const¶
Return number of vessel groups currently detected.
-
inline unsigned int get_n_coefficients() const¶
Return number of retained Fourier coefficients.
-
inline unsigned int get_offset_coefficients() const¶
Return coefficient offset used for harmonic indexing.
-
inline unsigned int get_inclusions_in_vessel(unsigned int vessel_id) const¶
Return number of inclusions mapped to a given vessel id.
-
inline std::vector<unsigned int> get_selected_coefficients() const¶
-
inline void set_n_q_points(unsigned int n_q)¶
-
inline void set_n_coefficients(unsigned int n_coefficients)¶
Override the number of Fourier coefficients per inclusion.
-
inline void set_fourier_coefficients(std::vector<unsigned int> temp)¶
-
inline void build_segment_index_vector()¶
-
inline unsigned int get_segment_index(const types::global_dof_index &quadrature_id) const¶
Public Members
-
double modulation_frequency = 0.0¶
Frequency used to modulate inclusion boundary data in time.
-
Particles::ParticleHandler<spacedim> inclusions_as_particles¶
Particle representation of inclusion quadrature points.
-
std::vector<std::vector<double>> inclusions¶
-
std::string data_file = ""¶
Optional ASCII file containing coefficient data.
-
std::vector<std::vector<double>> inclusions_data¶
Per-inclusion coefficient vectors used by coupling assembly.
-
std::vector<double> reference_inclusion_data¶
Default coefficients replicated when no data file is provided.
-
std::map<unsigned int, std::vector<types::global_dof_index>> map_vessel_inclusions¶
Mapping from vessel id to indices of associated inclusions.
-
bool cluster_with_segments = false¶
Switch for the clustering of inclusions dofs into segments.
Friends
- friend class CouplingOperator
-
inline Inclusions(const unsigned int n_vector_components = 1)¶