Reduced Lagrange Multipliers
 
Loading...
Searching...
No Matches
tensor_product_space.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2024 by Luca Heltai
4//
5// This file is part of the reduced_lagrange_multipliers application, based on
6// the deal.II library.
7//
8// The reduced_lagrange_multipliers application is free software; you can use
9// it, redistribute it, and/or modify it under the terms of the Apache-2.0
10// License WITH LLVM-exception as published by the Free Software Foundation;
11// either version 3.0 of the License, or (at your option) any later version. The
12// full text of the license can be found in the file LICENSE.md at the top level
13// of the reduced_lagrange_multipliers distribution.
14//
15// ---------------------------------------------------------------------
16
17#ifndef tensor_product_space_h
18#define tensor_product_space_h
19
20#include <deal.II/base/mpi.h>
23
26
28
31
33#include <deal.II/grid/tria.h>
34
36
39
41
42using namespace dealii;
43
56template <int reduced_dim, int dim, int spacedim, int n_components>
58{
66
73 static constexpr int cross_section_dim = dim - reduced_dim;
74
84
91 unsigned int fe_degree = 1;
92
101 unsigned int n_q_points = 0;
102
106 double thickness = 0.01;
107
114 std::string thickness_field_name = "";
115
119 std::string reduced_grid_name = "";
120};
121
122
146template <int reduced_dim, int dim, int spacedim, int n_components>
148{
149public:
159 &par,
160 MPI_Comm mpi_communicator = MPI_COMM_WORLD);
161
168 static constexpr int cross_section_dim = dim - reduced_dim;
169
173 std::function<void(Triangulation<reduced_dim, spacedim> &)>
175
185 std::function<void(
187 set_partitioner = [](auto &) {};
188
195 void
196 initialize();
197
203 const ReferenceCrossSection<dim - reduced_dim, spacedim, n_components> &
205
206
207 void
209
216 get_dof_handler() const;
217
218 const std::vector<Point<spacedim>> &
220
221 const std::vector<std::vector<double>> &
223
224 const std::vector<Point<spacedim>> &
226
227 const std::vector<std::vector<double>> &
229
230 const std::vector<std::vector<double>> &
232
239 void
241 const std::map<unsigned int, IndexSet> &remote_q_point_indices);
242
243
244 const std::vector<types::global_dof_index> &
246
247
256 std::tuple<unsigned int, unsigned int, unsigned int>
257 particle_id_to_cell_and_qpoint_indices(const unsigned int qpoint_index) const;
258
259
267 locally_owned_qpoints() const;
268
269
278
284 auto
285 get_quadrature() const -> const QGauss<reduced_dim> &;
286
287 void
289
290 const parallel::fullydistributed::Triangulation<reduced_dim, spacedim> &
291 get_triangulation() const;
292
293 double
294 get_scaling(const unsigned int) const;
295
296 const LinearAlgebra::distributed::Vector<double> &
297 get_properties() const;
298
299 const DoFHandler<reduced_dim, spacedim> &
300 get_properties_dh() const;
301
302 DoFHandler<reduced_dim, spacedim> &
304
305 const std::vector<std::string> &
306 get_properties_names() const;
307
308 std::vector<std::string> &
310
311protected:
318 void
319 setup_dofs();
320
325 std::map<unsigned int, IndexSet>
327 const std::map<unsigned int, IndexSet> &local_q_point_indices) const;
328
336
343 const TensorProductSpaceParameters<reduced_dim, dim, spacedim, n_components>
345
351 ReferenceCrossSection<cross_section_dim, spacedim, n_components>
353
359 parallel::fullydistributed::Triangulation<reduced_dim, spacedim>
361
368 FESystem<reduced_dim, spacedim> fe;
369
374
381 DoFHandler<reduced_dim, spacedim> dof_handler;
382
386 std::map<types::global_cell_index, std::vector<types::global_dof_index>>
388
389 std::vector<Point<spacedim>> all_qpoints;
390 std::vector<std::vector<double>> all_weights;
391
392 std::vector<Point<spacedim>> reduced_qpoints;
393 std::vector<std::vector<double>> reduced_weights;
394
398 LinearAlgebra::distributed::Vector<double> properties;
399
403 DoFHandler<reduced_dim, spacedim> properties_dh;
404
408 std::vector<std::string> properties_names;
409};
410
411
412// Template specializations for the TensorProductSpaceParameters
413
414
415
416#endif // tensor_product_space_h
ParameterAcceptor(const std::string &section_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
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
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
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
unsigned int cell_index
unsigned int global_cell_index
Parameter configuration for a ReferenceCrossSection.
ReferenceCrossSectionParameters< cross_section_dim, spacedim, n_components > section
std::string reduced_grid_name
Name of the grid to read from a file.
static constexpr int cross_section_dim