Template Class TensorProductSpace

Class Documentation

template<int reduced_dim, int dim, int spacedim, int n_components>
class TensorProductSpace

A class representing a tensor product space combining a lower-dimensional triangulation and a reference cross-section.

A class representing a tensor product space for reduced-dimensional problems.

This class provides functionality to work with tensor product spaces, including the initialization of reduced grids, handling degrees of freedom (DoFs), and managing reference cross-sections for reduced-dimensional domains.

Template Parameters:
  • reduced_dim – The dimension of the reduced triangulation.

  • dim – The dimension of the full-order object.

  • spacedim – The dimension of the ambient space

  • n_components – The number of components of the problem.

  • reduced_dim – The reduced dimension of the space.

  • dim – The full dimension of the space.

  • spacedim – The spatial dimension.

  • n_components – The number of components in the system.

Public Functions

TensorProductSpace(const TensorProductSpaceParameters<reduced_dim, dim, spacedim, n_components> &par, MPI_Comm mpi_communicator = MPI_COMM_WORLD)

Constructor for the TensorProductSpace class.

Initializes the tensor product space using the provided parameters.

Parameters:
  • par – The parameters defining the tensor product space.

  • mpi_communicator – Communicator used for distributed setup.

void initialize()

Initializes the tensor product space.

This function sets up the necessary components, such as the DoFHandler and finite element system.

const ReferenceCrossSection<dim - reduced_dim, spacedim, n_components> &get_reference_cross_section() const

Retrieves the reference cross-section.

Returns:

A constant reference to the reference cross-section object.

void make_reduced_grid_and_properties()

Build reduced grid and read optional reduced properties from file.

const DoFHandler<reduced_dim, spacedim> &get_dof_handler() const

Retrieves the DoFHandler for the reduced domain.

Returns:

A constant reference to the DoFHandler object.

const std::vector<Point<spacedim>> &get_locally_owned_qpoints() const

Return lifted quadrature points owned by this rank.

const std::vector<std::vector<double>> &get_locally_owned_weights() const

Return lifted quadrature weights owned by this rank.

const std::vector<Point<spacedim>> &get_locally_owned_reduced_qpoints() const

Return reduced-manifold quadrature points owned by this rank.

const std::vector<std::vector<double>> &get_locally_owned_reduced_weights() const

Return reduced-manifold quadrature weights owned by this rank.

const std::vector<std::vector<double>> &get_locally_owned_section_measure() const

Return reduced cross-section measures associated with local points.

void update_local_dof_indices(const std::map<unsigned int, IndexSet> &remote_q_point_indices)

Update the relevant local dof_indices.

After inserting global particles, this function updates the indices that are required to assemble the coupling matrix.

const std::vector<types::global_dof_index> &get_dof_indices(const types::global_cell_index cell_index) const

Return DoF indices associated with one reduced cell.

std::tuple<unsigned int, unsigned int, unsigned int> particle_id_to_cell_and_qpoint_indices(const unsigned int qpoint_index) const

Convert a global particle id to a global cell index, and the local quadrature indices on the reduce triangulation and on the cross-section.

Parameters:

qpoint_index – The global quadrature-point index.

Returns:

std::tuple<unsigned int, unsigned int, unsigned int> cell_index, q_index, qpoint_index_in_section

IndexSet locally_owned_qpoints() const

Return the indices of the quadrature points that are locally owned by the reduced domain.

Returns:

IndexSet

IndexSet locally_relevant_indices() const

Return the indices of the cells that are required to assemble the coupling matrix.

Returns:

IndexSet

auto get_quadrature() const -> const Quadrature<reduced_dim>&

Retrieve the quadrature formula used in the reduced domain.

Returns:

A constant reference to the quadrature formula.

void compute_points_and_weights()

Compute and cache lifted/reduced quadrature points and weights.

const parallel::fullydistributed::Triangulation<reduced_dim, spacedim> &get_triangulation() const

Return distributed reduced triangulation.

double get_scaling(const unsigned int) const

Return the scaling associated with one reduced cell.

const LinearAlgebra::distributed::Vector<double> &get_properties() const

Return interpolated reduced-grid properties.

const DoFHandler<reduced_dim, spacedim> &get_properties_dh() const

Return DoFHandler used to represent reduced-grid properties.

DoFHandler<reduced_dim, spacedim> &get_properties_dh()

Mutable access to the property DoFHandler.

const std::vector<std::string> &get_properties_names() const

Return names of reduced-grid property fields.

std::vector<std::string> &get_properties_names()

Mutable access to property field names.

Public Static Attributes

static constexpr int cross_section_dim = dim - reduced_dim

The dimension of the cross-section of the reduced domain.

This is computed as the difference between the full dimension and the reduced dimension.

Protected Functions

void setup_dofs()

Sets up the degrees of freedom (DoFs) for the reduced domain.

This function initializes the DoFHandler and associates it with the finite element system.

std::map<unsigned int, IndexSet> local_q_point_indices_to_global_cell_indices(const std::map<unsigned int, IndexSet> &local_q_point_indices) const

Given a map of processor to local quadrature point indices, return a map of processor to the corresponding global cell indices.

Protected Attributes

MPI_Comm mpi_communicator

The MPI communicator for parallel processing.

This object is used to manage communication between different processes in a parallel environment.

const TensorProductSpaceParameters<reduced_dim, dim, spacedim, n_components> &par

The parameters defining the tensor product space.

This object contains all the necessary configuration for the tensor product space.

ReferenceCrossSection<cross_section_dim, spacedim, n_components> reference_cross_section

The reference cross-section for the reduced domain.

This object represents the cross-section of the reduced-dimensional domain.

parallel::fullydistributed::Triangulation<reduced_dim, spacedim> triangulation

The triangulation representing the reduced domain.

This object holds the mesh for the reduced-dimensional domain.

FESystem<reduced_dim, spacedim> fe

The finite element system used for the reduced domain.

This object defines the finite element basis functions for the reduced-dimensional domain.

Quadrature<reduced_dim> quadrature_formula

The quadrature formula used for integration in the reduced domain.

DoFHandler<reduced_dim, spacedim> dof_handler

The DoFHandler for the reduced domain.

This object manages the degrees of freedom for the reduced-dimensional domain.

std::map<types::global_cell_index, std::vector<types::global_dof_index>> global_cell_to_dof_indices

Mapping from global cell index to dof indices.

std::vector<Point<spacedim>> all_qpoints

All quadrature points lifted in ambient coordinates.

std::vector<std::vector<double>> all_weights

Weights associated with all lifted quadrature points.

std::vector<Point<spacedim>> reduced_qpoints

Quadrature points on the reduced manifold.

std::vector<std::vector<double>> reduced_weights

Weights associated with reduced-manifold quadrature points.

LinearAlgebra::distributed::Vector<double> properties

The properties of the inclusion.

DoFHandler<reduced_dim, spacedim> properties_dh

The finite element system used for the properties of the inclusion.

std::vector<std::string> properties_names

The names of the properties stored in the input file.