32template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
45template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
55 ,
fe(
FE_Q<reduced_dim, spacedim>(
par.fe_degree),
64template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
84template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
92template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
101 VTKUtils::read_vtk(
par.reduced_grid_name,
102 serial_properties_dh,
106 deallog <<
"Read VTK file: " <<
par.reduced_grid_name
107 <<
", properties norm: " << serial_properties.
l2_norm() << std::endl;
121 <<
" has no locally owned cells." << std::endl;
131 VTKUtils::serial_vector_to_distributed_vector(serial_properties_dh,
139 deallog <<
"PROPERTIES ARE ZERO. DO NOT REFINE INPUT GRID" << std::endl;
146 const auto block_indices = VTKUtils::get_block_indices(properties_fe);
148 for (
unsigned int i = 0; i < block_indices.size(); ++i)
151 deallog <<
"Property name: " << name <<
", block index: " << i
152 <<
", block size: " << block_indices.block_size(i)
153 <<
", block start: " << block_indices.block_start(i) << std::endl;
156 deallog <<
"Serial properties norm: " << serial_properties.
l2_norm()
158 AssertDimension(block_indices.total_size(), properties_fe.n_components());
162template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
176template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
177const std::vector<Point<spacedim>> &
182 const int global_qpoints =
188 "No quadrature points exist across all MPI ranks. You must call compute_points_and_weights() first"));
192template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
193const std::vector<std::vector<double>> &
198 const int global_weights =
201 ExcMessage(
"No weights exist across all MPI ranks. You must call"
202 " compute_points_and_weights() first"));
206template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
207const std::vector<Point<spacedim>> &
212 const int global_reduced_qpoints =
215 global_reduced_qpoints > 0,
217 "No reduced quadrature points exist across all MPI ranks. You must call"
218 " compute_points_and_weights() first"));
222template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
223const std::vector<std::vector<double>> &
228 const int global_reduced_weights =
232 "No reduced weights exist across all MPI ranks. You must call"
233 " compute_points_and_weights() first"));
238template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
242 const std::map<unsigned int, IndexSet> &remote_q_point_indices)
244 auto global_cell_indices =
248 std::map<types::global_cell_index, std::vector<types::global_dof_index>>>
251 for (
const auto &[proc, cell_indices] : global_cell_indices)
252 for (
const auto &
id : cell_indices)
256 auto local_dof_indices =
259 for (
const auto &[proc, cell_indices] : local_dof_indices)
266template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
267const std::vector<types::global_dof_index> &
273 ExcMessage(
"Cell index not found in global cell to dof indices."));
279template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
280std::tuple<unsigned int, unsigned int, unsigned int>
291 const unsigned int qpoint_index_in_cell =
294 const unsigned int qpoint_index_in_section =
302 qpoint_index_in_cell,
303 qpoint_index_in_section);
307template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
308std::map<unsigned int, IndexSet>
311 const std::map<unsigned int, IndexSet> &remote_q_point_indices)
const
313 std::map<unsigned int, IndexSet> cell_indices;
317 ->locally_owned_range();
319 auto local_q_point_indices =
322 for (
const auto &[proc, qpoint_indices] : local_q_point_indices)
325 for (
const auto &qpoint_index : qpoint_indices)
329 cell_indices_for_proc.add_index(
332 cell_indices_for_proc.compress();
333 cell_indices[proc] = cell_indices_for_proc;
340template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
346 for (
const auto &cell :
dof_handler.active_cell_iterators())
347 if (cell->is_locally_owned())
349 std::vector<types::global_dof_index> dof_indices(
fe.n_dofs_per_cell());
350 cell->get_dof_indices(dof_indices);
356template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
364 ->locally_owned_range();
369 const unsigned int n_qpoints_per_cell =
375 return locally_owned_qpoints_set;
378template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
385 ->locally_owned_range();
392template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
395 const -> const
QGauss<reduced_dim> &
400template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
432 const unsigned int thickness_field_index =
436 par.thickness_field_name));
438 const auto thickness_start =
441 VTKUtils::get_block_indices(properties_fe)
442 .block_start(thickness_field_index);
449 for (
const auto &cell :
triangulation.active_cell_iterators())
450 if (cell->is_locally_owned())
459 properties_fe_values.
reinit(
460 cell->as_dof_handler_iterator(this->properties_dh));
468 for (
const auto &
w : JxW)
472 if constexpr (reduced_dim == 1)
473 new_vertical = cell->vertex(1) - cell->vertex(0);
477 const auto &qpoint = qpoints[q];
478 if constexpr (reduced_dim == 2)
481 auto cross_section_qpoints =
483 qpoint, new_vertical, thickness_values[q]);
486 cross_section_qpoints.get_points().begin(),
487 cross_section_qpoints.get_points().end());
489 for (
const auto &
w : cross_section_qpoints.get_weights())
495template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
503template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
506 const unsigned int)
const
508 return std::pow(
par.thickness, -((dim - reduced_dim) / 2.0));
511template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
519template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
527template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
535template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
536const std::vector<std::string> &
543template <
int reduced_dim,
int dim,
int spacedim,
int n_components>
544std::vector<std::string> &
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
types::global_dof_index n_dofs() const
const std::vector< double > & get_JxW_values() const
void get_function_values(const ReadVector< Number > &fe_function, std::vector< Number > &values) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const Tensor< 1, spacedim > & normal_vector(const unsigned int q_point) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
double JxW(const unsigned int q_point) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
IndexSet tensor_product(const IndexSet &other) const
void add_index(const size_type index)
size_type nth_index_in_set(const size_type local_index) const
ParameterAcceptor(const std::string §ion_name="")
void add_parameter(const std::string &entry, ParameterType ¶meter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern())
real_type l2_norm() const
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
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
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::function< void(Triangulation< reduced_dim, spacedim > &)> preprocess_serial_triangulation
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
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define AssertThrow(cond, exc)
IndexSet complete_index_set(const IndexSet::size_type N)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T sum(const T &t, const MPI_Comm mpi_communicator)
std::map< unsigned int, T > some_to_some(const MPI_Comm comm, const std::map< unsigned int, T > &objects_to_send)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
constexpr unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
unsigned int global_cell_index
std::string thickness_field_name
TensorProductSpaceParameters()
std::string reduced_grid_name
Name of the grid to read from a file.