17#ifndef tensor_product_space_h
18#define tensor_product_space_h
56template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
146template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
212 const std::vector<Point<spacedim>> &
215 const std::vector<std::vector<double>> &
218 const std::vector<Point<spacedim>> &
221 const std::vector<std::vector<double>> &
224 const std::vector<std::vector<double>> &
235 const std::map<unsigned int, IndexSet> &remote_q_point_indices);
238 const std::vector<types::global_dof_index> &
250 std::tuple<unsigned int, unsigned int, unsigned int>
299 const
std::vector<
std::
string> &
302 std::vector<
std::
string> &
321 const
std::map<
unsigned int,
IndexSet> &local_q_point_indices) const;
ParameterAcceptor(const std::string §ion_name="")
Handles the construction and management of a reference inclusion geometry and its basis.
parallel::fullydistributed::Triangulation< reduced_dim, spacedim > triangulation
std::vector< Point< spacedim > > reduced_qpoints
std::vector< std::vector< double > > reduced_weights
MPI_Comm mpi_communicator
const LinearAlgebra::distributed::Vector< double > & get_properties() const
const parallel::fullydistributed::Triangulation< reduced_dim, spacedim > & get_triangulation() const
std::vector< std::vector< double > > all_weights
const std::vector< Point< spacedim > > & get_locally_owned_qpoints() const
TensorProductSpace(const TensorProductSpaceParameters< reduced_dim, dim, spacedim, n_components > &par, MPI_Comm mpi_communicator=MPI_COMM_WORLD)
DoFHandler< reduced_dim, spacedim > properties_dh
double get_scaling(const unsigned int) const
void compute_points_and_weights()
const ReferenceCrossSection< dim - reduced_dim, spacedim, n_components > & get_reference_cross_section() const
auto get_quadrature() const -> const QGauss< reduced_dim > &
std::vector< Point< spacedim > > all_qpoints
std::function< void(parallel::fullydistributed::Triangulation< reduced_dim, spacedim > &)> set_partitioner
const std::vector< std::vector< double > > & get_locally_owned_section_measure() const
std::map< types::global_cell_index, std::vector< types::global_dof_index > > global_cell_to_dof_indices
const DoFHandler< reduced_dim, spacedim > & get_dof_handler() const
void update_local_dof_indices(const std::map< unsigned int, IndexSet > &remote_q_point_indices)
const std::vector< types::global_dof_index > & get_dof_indices(const types::global_cell_index cell_index) const
static constexpr int cross_section_dim
const DoFHandler< reduced_dim, spacedim > & get_properties_dh() const
const std::vector< std::vector< double > > & get_locally_owned_reduced_weights() const
void make_reduced_grid_and_properties()
IndexSet locally_relevant_indices() const
std::vector< std::string > properties_names
IndexSet locally_owned_qpoints() const
std::tuple< unsigned int, unsigned int, unsigned int > particle_id_to_cell_and_qpoint_indices(const unsigned int qpoint_index) const
DoFHandler< reduced_dim, spacedim > dof_handler
const std::vector< std::string > & get_properties_names() const
const std::vector< Point< spacedim > > & get_locally_owned_reduced_qpoints() const
LinearAlgebra::distributed::Vector< double > properties
QGauss< reduced_dim > quadrature_formula
const TensorProductSpaceParameters< reduced_dim, dim, spacedim, n_components > & par
std::map< unsigned int, IndexSet > local_q_point_indices_to_global_cell_indices(const std::map< unsigned int, IndexSet > &local_q_point_indices) const
const std::vector< std::vector< double > > & get_locally_owned_weights() const
ReferenceCrossSection< cross_section_dim, spacedim, n_components > reference_cross_section
FESystem< reduced_dim, spacedim > fe
unsigned int global_cell_index
Parameter configuration for a ReferenceCrossSection.
std::string thickness_field_name
ReferenceCrossSectionParameters< cross_section_dim, spacedim, n_components > section
TensorProductSpaceParameters()
std::string reduced_grid_name
Name of the grid to read from a file.
static constexpr int cross_section_dim