/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FEM/assembler/cg.hh Source File#
|
DiFfRG
|
cg.hh
Go to the documentation of this file.
25 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "fe_derivatives", "fe_hessians">>(
90 solution_grad_interface[0].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
91 solution_grad_interface[1].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
92 solution_hess_interface[0].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
93 solution_hess_interface[1].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
164 template <typename Discretization_, typename Model_> class Assembler : public FEMAssembler<Discretization_, Model_>
186 virtual void reinit_vector(VectorType &vec) const override { vec.reinit(dof_handler.n_dofs()); }
201 MatrixCreator::create_mass_matrix(dof_handler, quadrature, mass_matrix, (Function<dim, NumberType> *)nullptr,
251 virtual const SparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
260 virtual void refinement_indicator(Vector<double> &indicator, const VectorType &solution_global) override
266 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
297 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
301 virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot,
309 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
334 copy_data.cell_residual(i) += weight * JxW[q_index] * fe_v.shape_value_component(i, q_index, comp[i]) *
347 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
351 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
366 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
396 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
399 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
409 copy_data.cell_residual(i) += weight_mass * JxW[q_index] * fe_v.shape_value_component(i, q_index, ci) *
414 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
441 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
458 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces;
461 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
466 virtual void jacobian_mass(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
474 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
520 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
525 virtual void jacobian(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global, NumberType weight,
542 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
567 SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()> j_grad_flux;
568 SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()> j_hess_flux;
569 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
574 SimpleMatrix<NumberType, Components::count_fe_functions(), Components::count_extractors()> j_extr_source;
582 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
585 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
588 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
592 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
620 contribution += JxW[q_index] * sv[j] * sv_i * (alpha * j_mass_dot(ci, cj) + beta * j_mass(ci, cj));
630 tbb::blocked_range2d<uint>(0, n_dofs, 0, n_dofs), [&](const tbb::blocked_range2d<uint> &range) {
649 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
676 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
682 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
685 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
688 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
692 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
726 tbb::blocked_range2d<uint>(0, n_dofs, 0, n_dofs), [&](const tbb::blocked_range2d<uint> &range) {
727 do_bnd_work(range.rows().begin(), range.rows().end(), range.cols().begin(), range.cols().end());
747 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices.size(), extractor_dof_indices.size());
749 constraints.distribute_local_to_global(extractor_dependence, c.local_dof_indices, extractor_dof_indices,
756 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces;
759 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
768 ss << " Reinit: " << average_time_reinit() * 1000 << "ms (" << num_reinits() << ")" << std::endl;
769 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
771 ss << " Jacobian: " << average_time_jacobian_assembly() * 1000 << "ms (" << num_jacobians() << ")"
The basic assembler that can be used for any standard CG scheme with flux and source.
Definition cg.hh:165
double average_time_residual_assembly()
Definition cg.hh:786
double average_time_jacobian_assembly()
Definition cg.hh:796
virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight) override
Definition cg.hh:301
virtual void refinement_indicator(Vector< double > &indicator, const VectorType &solution_global) override
refinement indicator for adaptivity. Only calls the model's cell_indicator function,...
Definition cg.hh:260
virtual void rebuild_jacobian_sparsity() override
Definition cg.hh:235
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition cg.hh:188
SparsityPattern sparsity_pattern_mass
Definition cg.hh:818
Assembler(Discretization &discretization, Model &model, const JSONValue &json)
Definition cg.hh:177
virtual const SparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition cg.hh:247
virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType weight_mass, const VectorType &variables=VectorType()) override
Definition cg.hh:351
SparsityPattern sparsity_pattern_jacobian
Definition cg.hh:819
virtual const SparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition cg.hh:251
virtual void jacobian(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType alpha, NumberType beta, const VectorType &variables=VectorType()) override
Definition cg.hh:525
std::vector< types::global_dof_index > extractor_dof_indices
Definition common.hh:355
virtual void jacobian_mass(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha, NumberType beta) override
Definition cg.hh:466
virtual void reinit_vector(VectorType &vec) const override
Definition cg.hh:186
The basic assembler that can be used for any standard CG scheme with flux and source.
Definition common.hh:39
FullMatrix< NumberType > extractor_jacobian
Definition common.hh:351
void extract(std::array< NumberType, Components::count_extractors()> &data, const VectorType &solution_global, const VectorType &variables, bool search_EoM, bool set_EoM, bool postprocess) const
Definition common.hh:200
bool jacobian_extractors(FullMatrix< NumberType > &extractor_jacobian, const VectorType &solution_global, const VectorType &variables)
Definition common.hh:237
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition common.hh:108
std::vector< types::global_dof_index > extractor_dof_indices
Definition common.hh:355
A simple NxM-matrix class, which is used for cell-wise Jacobians.
Definition tuples.hh:156
Definition complex_math.hh:10
Definition cg.hh:149
std::array< double, 2 > values
Definition cg.hh:151
std::array< uint, 2 > cell_indices
Definition cg.hh:150
Definition cg.hh:148
std::vector< CopyFaceData_I > face_data
Definition cg.hh:153
Definition cg.hh:134
FullMatrix< NumberType > extractor_cell_jacobian
Definition cg.hh:136
std::vector< types::global_dof_index > local_dof_indices
Definition cg.hh:137
FullMatrix< NumberType > cell_jacobian
Definition cg.hh:135
void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
Definition cg.hh:139
Definition cg.hh:122
void reinit(const Iterator &cell, uint dofs_per_cell)
Definition cg.hh:126
Vector< NumberType > cell_residual
Definition cg.hh:123
std::vector< types::global_dof_index > local_dof_indices
Definition cg.hh:124
Class to hold data for each assembly thread, i.e. FEValues for cells, interfaces, as well as pre-allo...
Definition cg.hh:37
std::vector< Vector< NumberType > > solution
Definition cg.hh:106
array< std::vector< std::vector< Tensor< 1, dim, NumberType > > >, 2 > solution_grad_interface
Definition cg.hh:111
array< std::vector< std::vector< Tensor< 2, dim, NumberType > > >, 2 > solution_hess_interface
Definition cg.hh:112
ScratchData(const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const dealii::Quadrature< dim > &quadrature, const dealii::Quadrature< dim - 1 > &quadrature_face, const UpdateFlags update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values|update_hessians, const UpdateFlags interface_update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values|update_normal_vectors|update_hessians)
Definition cg.hh:41
std::vector< Tensor< 1, dim > > cached_shape_grads
Definition cg.hh:118
std::vector< std::vector< Tensor< 1, dim, NumberType > > > solution_grad
Definition cg.hh:107
std::vector< double > cached_shape_values
Definition cg.hh:117
FEInterfaceValues< dim > fe_interface_values
Definition cg.hh:104
ScratchData(const ScratchData< Discretization > &scratch_data)
Definition cg.hh:72
std::vector< Tensor< 2, dim > > cached_shape_hessians
Definition cg.hh:119
array< std::vector< Vector< NumberType > >, 2 > solution_interface
Definition cg.hh:110
std::vector< std::vector< Tensor< 2, dim, NumberType > > > solution_hess
Definition cg.hh:108
std::vector< Vector< NumberType > > solution_dot
Definition cg.hh:109
typename Discretization::NumberType NumberType
Definition cg.hh:39
Definition tuples.hh:34
A class to store a tuple with elements that can be accessed by name. The names are stored as FixedStr...
Definition tuples.hh:56
Generated by