PRISMS-PF Manual v3.0-pre
Loading...
Searching...
No Matches
FloodFiller< dim, degree > Class Template Reference

#include <FloodFiller.h>

Public Member Functions

 FloodFiller (dealii::FESystem< dim > &_fe, dealii::QGaussLobatto< dim > _quadrature)
 
void calcGrainSets (dealii::FESystem< dim > &finite_element, dealii::DoFHandler< dim > &dof_handler, dealii::LinearAlgebra::distributed::Vector< double > *solution_field, double threshold_lower, double threshold_upper, int min_id, unsigned int order_parameter_index, std::vector< GrainSet< dim > > &grain_sets)
 

Protected Member Functions

template<typename T >
void recursiveFloodFill (T cell, T cell_end, dealii::LinearAlgebra::distributed::Vector< double > *solution_field, double threshold_lower, double threshold_upper, int min_id, unsigned int &grain_index, std::vector< GrainSet< dim > > &grain_sets, bool &grain_assigned)
 
void createGlobalGrainSetList (std::vector< GrainSet< dim > > &grain_sets) const
 
void mergeSplitGrains (std::vector< GrainSet< dim > > &grain_sets) const
 

Protected Attributes

dealii::QGaussLobatto< dim > quadrature
 
const unsigned int num_quad_points
 
const unsigned int dofs_per_cell
 
dealii::FESystem< dim > * fe
 

Detailed Description

template<int dim, int degree>
class FloodFiller< dim, degree >

This class uses a recursive flood filling algorithm to find connected bodies in a field, given a threshold. The MPI communication methods are similar to those in parallelNucleationList.

Constructor & Destructor Documentation

◆ FloodFiller()

template<int dim, int degree>
FloodFiller< dim, degree >::FloodFiller ( dealii::FESystem< dim > &  _fe,
dealii::QGaussLobatto< dim >  _quadrature 
)
inline

Constructor.

Member Function Documentation

◆ calcGrainSets()

template<int dim, int degree>
void FloodFiller< dim, degree >::calcGrainSets ( dealii::FESystem< dim > &  finite_element,
dealii::DoFHandler< dim > &  dof_handler,
dealii::LinearAlgebra::distributed::Vector< double > *  solution_field,
double  threshold_lower,
double  threshold_upper,
int  min_id,
unsigned int  order_parameter_index,
std::vector< GrainSet< dim > > &  grain_sets 
)

The primary external interface. This method takes in information about the mesh/field and outputs a vector of GrainSet objects.

◆ createGlobalGrainSetList()

template<int dim, int degree>
void FloodFiller< dim, degree >::createGlobalGrainSetList ( std::vector< GrainSet< dim > > &  grain_sets) const
protected

The method to merge the grain sets from all the processors.

◆ mergeSplitGrains()

template<int dim, int degree>
void FloodFiller< dim, degree >::mergeSplitGrains ( std::vector< GrainSet< dim > > &  grain_sets) const
protected

Checks to see if grains found on different processors are parts of a larger grain. If so, it merges the grain_sets entries.

◆ recursiveFloodFill()

template<int dim, int degree>
template<typename T >
void FloodFiller< dim, degree >::recursiveFloodFill ( cell,
cell_end,
dealii::LinearAlgebra::distributed::Vector< double > *  solution_field,
double  threshold_lower,
double  threshold_upper,
int  min_id,
unsigned int &  grain_index,
std::vector< GrainSet< dim > > &  grain_sets,
bool &  grain_assigned 
)
protected

The actual recursive flood fill method.

Member Data Documentation

◆ dofs_per_cell

template<int dim, int degree>
const unsigned int FloodFiller< dim, degree >::dofs_per_cell
protected

The number of degrees of freedom per cell.

◆ fe

template<int dim, int degree>
dealii::FESystem<dim>* FloodFiller< dim, degree >::fe
protected

The deal.II finite element object, set in the constructor.

◆ num_quad_points

template<int dim, int degree>
const unsigned int FloodFiller< dim, degree >::num_quad_points
protected

The number of quadrature points per cell.

◆ quadrature

template<int dim, int degree>
dealii::QGaussLobatto<dim> FloodFiller< dim, degree >::quadrature
protected

The quadrature used to calculate the element-wise value of the solution field.


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