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

This is the abstract base class for the matrix free implementation of parabolic and elliptic BVP's, and supports MPI, threads and vectorization (Hybrid Parallel). This class contains the parallel data structures, mesh (referred to as triangulation), parallel degrees of freedom distribution, constraints, and general utility methods. More...

#include <matrixFreePDE.h>

Inheritance diagram for MatrixFreePDE< dim, degree >:

Public Member Functions

 MatrixFreePDE (userInputParameters< dim >)
 Class contructor.
 
 ~MatrixFreePDE () override
 Class destructor.
 
virtual void init ()
 Initializes the mesh, degrees of freedom, constraints and data structures using the user provided inputs in the application parameters file.
 
virtual void create_triangulation (parallel::distributed::Triangulation< dim > &tria) const
 Create the mesh with the user provided domain sizes.
 
void initForTests (std::vector< Field< dim > > _fields)
 Initializes the data structures for enabling unit tests.
 
void solve ()
 This method implements the time stepping algorithm and invokes the solveIncrement() method.
 
void vmult (dealii::LinearAlgebra::distributed::Vector< double > &dst, const dealii::LinearAlgebra::distributed::Vector< double > &src) const
 This method essentially converts the MatrixFreePDE object into a matrix object which can be used with matrix free iterative solvers. Provides the A*x functionality for solving the system of equations Ax=b.
 
void buildFields ()
 Create the vector of all physical fields in the problem.
 
virtual void setInitialCondition (const Point< dim > &p, const unsigned int index, double &scalar_IC, Vector< double > &vector_IC)=0
 Set the initial condition for all fields. This function is overriden in each application.
 
virtual void setNonUniformDirichletBCs (const Point< dim > &p, const unsigned int index, const unsigned int direction, const double time, double &scalar_BC, Vector< double > &vector_BC)=0
 Set the spatially or temporally non-uniform boundary conditions. This function is overriden in each application.
 

Public Attributes

std::vector< Field< dim > > fields
 Vector of all the physical fields in the problem. Fields are identified by dimentionality (SCALAR/VECTOR), the kind of PDE (ELLIPTIC/PARABOLIC) used to compute them and a character identifier (e.g.: "c" for composition) which is used to write the fields to the output files.
 
ConditionalOStream pcout
 Parallel message stream.
 

Protected Member Functions

void reinit ()
 
void reassignGrains ()
 
virtual void solveIncrement (bool skip_time_dependent)
 
void outputResults ()
 
void computeInvM ()
 
void updateExplicitSolution (unsigned int fieldIndex)
 
bool updateImplicitSolution (unsigned int fieldIndex, unsigned int nonlinear_iteration_index)
 
void applyBCs (unsigned int fieldIndex)
 
void compute_element_volume ()
 Compute element volume for the triangulation.
 
void computeExplicitRHS ()
 
void computeNonexplicitRHS ()
 
void getLHS (const MatrixFree< dim, double > &data, dealii::LinearAlgebra::distributed::Vector< double > &dst, const dealii::LinearAlgebra::distributed::Vector< double > &src, const std::pair< unsigned int, unsigned int > &cell_range) const
 
void getLaplaceLHS (const MatrixFree< dim, double > &data, dealii::LinearAlgebra::distributed::Vector< double > &dst, const dealii::LinearAlgebra::distributed::Vector< double > &src, const std::pair< unsigned int, unsigned int > &cell_range) const
 
void setNonlinearEqInitialGuess ()
 
void computeLaplaceRHS (unsigned int fieldIndex)
 
void getLaplaceRHS (const MatrixFree< dim, double > &data, dealii::LinearAlgebra::distributed::Vector< double > &dst, const dealii::LinearAlgebra::distributed::Vector< double > &src, const std::pair< unsigned int, unsigned int > &cell_range) const
 
void getExplicitRHS (const MatrixFree< dim, double > &data, std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > &dst, const std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > &src, const std::pair< unsigned int, unsigned int > &cell_range) const
 
void getNonexplicitRHS (const MatrixFree< dim, double > &data, std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > &dst, const std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > &src, const std::pair< unsigned int, unsigned int > &cell_range) const
 
virtual void explicitEquationRHS (variableContainer< dim, degree, VectorizedArray< double > > &variable_list, const Point< dim, VectorizedArray< double > > q_point_loc, const VectorizedArray< double > element_volume) const =0
 
virtual void nonExplicitEquationRHS (variableContainer< dim, degree, VectorizedArray< double > > &variable_list, const Point< dim, VectorizedArray< double > > q_point_loc, const VectorizedArray< double > element_volume) const =0
 
virtual void equationLHS (variableContainer< dim, degree, VectorizedArray< double > > &variable_list, const Point< dim, VectorizedArray< double > > q_point_loc, const VectorizedArray< double > element_volume) const =0
 
virtual void postProcessedFields (const variableContainer< dim, degree, VectorizedArray< double > > &variable_list, variableContainer< dim, degree, VectorizedArray< double > > &pp_variable_list, const Point< dim, VectorizedArray< double > > q_point_loc, const VectorizedArray< double > element_volume) const
 
void computePostProcessedFields (std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > &postProcessedSet)
 
void getPostProcessedFields (const MatrixFree< dim, double > &data, std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > &dst, const std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > &src, const std::pair< unsigned int, unsigned int > &cell_range)
 
void markBoundaries (parallel::distributed::Triangulation< dim > &) const
 
void applyDirichletBCs ()
 
void applyNeumannBCs ()
 
void setPeriodicity ()
 
void setPeriodicityConstraints (AffineConstraints< double > *, const DoFHandler< dim > *) const
 
void set_rigid_body_mode_constraints (AffineConstraints< double > *constraints, const DoFHandler< dim > *dof_handler, const Point< dim > target_point=Point< dim >()) const
 Set constraints to pin the solution to 0 at a certain vertex. This is automatically done at the origin if no value terms are detected in your dependencies in a time_independent or implicit solve.
 
void applyInitialConditions ()
 
void save_checkpoint ()
 
void load_checkpoint_triangulation ()
 
void load_checkpoint_fields ()
 
void load_checkpoint_time_info ()
 
void move_file (const std::string &, const std::string &)
 
void verify_checkpoint_file_exists (const std::string &filename)
 
void updateNucleiList ()
 
std::vector< nucleus< dim > > getNewNuclei ()
 
void getLocalNucleiList (std::vector< nucleus< dim > > &newnuclei) const
 
void safetyCheckNewNuclei (std::vector< nucleus< dim > > newnuclei, std::vector< unsigned int > &conflict_ids)
 
void refineMeshNearNuclei (std::vector< nucleus< dim > > newnuclei)
 
double weightedDistanceFromNucleusCenter (const Point< dim, double > center, const std::vector< double > &semiaxes, const Point< dim, double > q_point_loc, const unsigned int var_index) const
 
VectorizedArray< double > weightedDistanceFromNucleusCenter (const Point< dim, double > center, const std::vector< double > &semiaxes, const Point< dim, VectorizedArray< double > > q_point_loc, const unsigned int var_index) const
 
virtual double getNucleationProbability (variableValueContainer, double, Point< dim >, unsigned int variable_index) const
 
unsigned int getFieldIndex (std::string _name)
 
void computeIntegral (double &integratedField, int index, std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > variableSet)
 

Protected Attributes

userInputParameters< dim > userInputs
 
unsigned int totalDOFs
 
const AttributesList & var_attributes
 
const AttributesList & pp_attributes
 
std::vector< SimplifiedGrainRepresentation< dim > > simplified_grain_representations
 
parallel::distributed::Triangulation< dim > triangulation
 
std::vector< FESystem< dim > * > FESet
 
std::vector< const AffineConstraints< double > * > constraintsDirichletSet
 
std::vector< const AffineConstraints< double > * > constraintsOtherSet
 
std::vector< const DoFHandler< dim > * > dofHandlersSet
 
std::vector< const IndexSet * > locally_relevant_dofsSet
 
std::vector< AffineConstraints< double > * > constraintsDirichletSet_nonconst
 
std::vector< AffineConstraints< double > * > constraintsOtherSet_nonconst
 
std::vector< DoFHandler< dim > * > dofHandlersSet_nonconst
 
std::vector< IndexSet * > locally_relevant_dofsSet_nonconst
 
std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > solutionSet
 
std::vector< dealii::LinearAlgebra::distributed::Vector< double > * > residualSet
 
std::vector< parallel::distributed::SolutionTransfer< dim, dealii::LinearAlgebra::distributed::Vector< double > > * > soltransSet
 
MatrixFree< dim, double > matrixFreeObject
 
dealii::LinearAlgebra::distributed::Vector< double > invMscalar
 
dealii::LinearAlgebra::distributed::Vector< double > invMvector
 
dealii::LinearAlgebra::distributed::Vector< double > dU_vector
 
dealii::LinearAlgebra::distributed::Vector< double > dU_scalar
 
unsigned int currentFieldIndex
 
dealii::AlignedVector< dealii::VectorizedArray< double > > element_volume
 Vector that stores element volumes.
 
bool generatingInitialGuess
 
std::vector< std::map< types::global_dof_index, double > * > valuesDirichletSet
 
std::vector< nucleus< dim > > nuclei
 
bool isTimeDependentBVP
 
bool isEllipticBVP
 
bool hasExplicitEquation
 
bool hasNonExplicitEquation
 
double currentTime
 
unsigned int currentIncrement
 
unsigned int currentOutput
 
unsigned int currentCheckpoint
 
unsigned int current_grain_reassignment
 
TimerOutput computing_timer
 
bool first_integrated_var_output_complete
 
double integrated_var
 
unsigned int integral_index
 
std::mutex assembler_lock
 
AdaptiveRefinement< dim, degree > AMR
 

Static Protected Attributes

static const unsigned int CIJ_tensor_size = 2 * dim - 1 + dim / 3
 

Detailed Description

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

This is the abstract base class for the matrix free implementation of parabolic and elliptic BVP's, and supports MPI, threads and vectorization (Hybrid Parallel). This class contains the parallel data structures, mesh (referred to as triangulation), parallel degrees of freedom distribution, constraints, and general utility methods.

All the physical models in this package inherit this base class.

Template Parameters
dimThe number of dimensions in the problem.
degreeThe polynomial degree of the shape functions.

Member Function Documentation

◆ applyDirichletBCs()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::applyDirichletBCs ( )
protected

Method for applying Dirichlet boundary conditions.

◆ applyNeumannBCs()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::applyNeumannBCs ( )
protected

Method for applying Neumann boundary conditions.

◆ create_triangulation()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::create_triangulation ( parallel::distributed::Triangulation< dim > &  tria) const
virtual

Create the mesh with the user provided domain sizes.

Parameters
triaTriangulation object. It must be empty prior to calling this function.

◆ initForTests()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::initForTests ( std::vector< Field< dim > >  _fields)

Initializes the data structures for enabling unit tests.

This method initializes the MatrixFreePDE object with a fixed geometry, discretization and other custom selected options specifically to help with unit tests, and should not be called in any of the physical models.

Parameters
_fieldsVector of PDE descriptions (e.g., scalar/vector) for each field.

◆ reassignGrains()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::reassignGrains ( )
protected

Method to reassign grains when multiple grains are stored in a single order parameter.

◆ set_rigid_body_mode_constraints()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::set_rigid_body_mode_constraints ( AffineConstraints< double > *  constraints,
const DoFHandler< dim > *  dof_handler,
const Point< dim >  target_point = Point<dim>() 
) const
protected

Set constraints to pin the solution to 0 at a certain vertex. This is automatically done at the origin if no value terms are detected in your dependencies in a time_independent or implicit solve.

Parameters
constraintsThe constraint set.
dof_handlerThe list of the degrees of freedom.
target_pointThe point where the solution is constrained. This is the origin by default.

◆ setInitialCondition()

template<int dim, int degree>
virtual void MatrixFreePDE< dim, degree >::setInitialCondition ( const Point< dim > &  p,
const unsigned int  index,
double &  scalar_IC,
Vector< double > &  vector_IC 
)
pure virtual

Set the initial condition for all fields. This function is overriden in each application.

Parameters
pThe point at which the initial condition is evaluated.
indexThe index of the field being evaluated.
scalar_ICReturn variable for scalar fields.
vector_ICReturn variable for vector fields.

◆ setNonUniformDirichletBCs()

template<int dim, int degree>
virtual void MatrixFreePDE< dim, degree >::setNonUniformDirichletBCs ( const Point< dim > &  p,
const unsigned int  index,
const unsigned int  direction,
const double  time,
double &  scalar_BC,
Vector< double > &  vector_BC 
)
pure virtual

Set the spatially or temporally non-uniform boundary conditions. This function is overriden in each application.

Parameters
pThe point at which the boundary condition is evaluated.
indexThe index of the field being evaluated.
directionThe face of the boundary condition.
timeThe time at which the boundary condition is evaluated.
scalar_BCReturn variable for scalar fields.
vector_BCReturn variable for vector fields.

◆ solveIncrement()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::solveIncrement ( bool  skip_time_dependent)
protectedvirtual

Method to solve each time increment of a time-dependent problem. For time-independent problems this method is called only once. This method solves for all the fields in a staggered manner (one after another) and also invokes the corresponding solvers: Explicit solver for Parabolic problems, Implicit (matrix-free) solver for Elliptic problems.

◆ vmult()

template<int dim, int degree>
void MatrixFreePDE< dim, degree >::vmult ( dealii::LinearAlgebra::distributed::Vector< double > &  dst,
const dealii::LinearAlgebra::distributed::Vector< double > &  src 
) const

This method essentially converts the MatrixFreePDE object into a matrix object which can be used with matrix free iterative solvers. Provides the A*x functionality for solving the system of equations Ax=b.

Parameters
dstThe destination vector.
srcThe source vector.

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