/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FEM/assembler/ddg.hh Source File#
|
DiFfRG
|
ddg.hh
Go to the documentation of this file.
25 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "fe_derivatives", "fe_hessians">>(
41 ScratchData(const Mapping<dim> &mapping, const FiniteElement<dim> &fe, const Quadrature<dim> &quadrature,
89 solution_grad_interface[0].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
90 solution_grad_interface[1].resize(q_face_size, std::vector<Tensor<1, dim, NumberType>>(n_components));
91 solution_hess_interface[0].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
92 solution_hess_interface[1].resize(q_face_size, std::vector<Tensor<2, dim, NumberType>>(n_components));
191 template <typename Discretization_, typename Model_> class Assembler : public FEMAssembler<Discretization_, Model_>
213 virtual void reinit_vector(VectorType &vec) const override { vec.reinit(dof_handler.n_dofs()); }
228 MatrixCreator::create_mass_matrix(dof_handler, quadrature, mass_matrix, (Function<dim, NumberType> *)nullptr,
258 virtual const SparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
266 virtual void refinement_indicator(Vector<double> &indicator, const VectorType &solution_global) override
272 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
297 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
349 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_own_interior_faces_once;
351 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
355 virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot,
363 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
403 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
407 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
422 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
452 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
455 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
471 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
498 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
505 scalar_product(numflux[component_i], normals[q_index])); // phi_i(x_q) * numflux(x_q, u_q) * n(x_q)
509 const auto face_worker = [&](const Iterator &cell, const uint &f, const unsigned int &sf, const Iterator &ncell,
544 comp[i] = cd_i[0] == numbers::invalid_unsigned_int ? fe.system_to_component_index(cd_i[1]).first
575 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
579 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
584 virtual void jacobian_mass(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
592 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
640 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
645 virtual void jacobian(SparseMatrix<NumberType> &jacobian, const VectorType &solution_global, NumberType weight,
663 const auto cell_worker = [&](const Iterator &cell, Scratch &scratch_data, CopyData ©_data) {
692 SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()> j_grad_flux;
693 SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()> j_hess_flux;
694 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
699 SimpleMatrix<NumberType, Components::count_fe_functions(), Components::count_extractors()> j_extr_source;
707 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
710 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
713 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
717 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
769 const auto boundary_worker = [&](const Iterator &cell, const uint &face_no, Scratch &scratch_data,
800 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
806 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
809 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
812 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
816 fe_tie(solution[q_index], solution_grad[q_index], solution_hess[q_index], extracted_data, variables));
861 const auto face_worker = [&](const Iterator &cell, const uint &f, const uint &sf, const Iterator &ncell,
911 array<SimpleMatrix<Tensor<1, dim, Tensor<1, dim, NumberType>>, Components::count_fe_functions()>, 2>
913 array<SimpleMatrix<Tensor<1, dim, Tensor<2, dim, NumberType>>, Components::count_fe_functions()>, 2>
947 fe_iv.get_fe_face_values(face_no_arr[k]).shape_value_component(local_dof_arr[k], q_index, comp[k]);
949 fe_iv.get_fe_face_values(face_no_arr[k]).shape_grad_component(local_dof_arr[k], q_index, comp[k]);
951 fe_iv.get_fe_face_values(face_no_arr[k]).shape_hessian_component(local_dof_arr[k], q_index, comp[k]);
996 FullMatrix<NumberType> extractor_dependence(cdf.joint_dof_indices.size(), extractor_dof_indices.size());
998 constraints.distribute_local_to_global(extractor_dependence, cdf.joint_dof_indices, extractor_dof_indices,
1003 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices.size(), extractor_dof_indices.size());
1005 constraints.distribute_local_to_global(extractor_dependence, c.local_dof_indices, extractor_dof_indices,
1012 MeshWorker::AssembleFlags flags = MeshWorker::assemble_own_cells | MeshWorker::assemble_boundary_faces |
1016 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
1025 ss << " Reinit: " << average_time_reinit() * 1000 << "ms (" << num_reinits() << ")" << std::endl;
1026 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
1028 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 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
Definition quadrature.hh:56
size_t size() const
A simple NxM-matrix class, which is used for cell-wise Jacobians.
Definition tuples.hh:156
The basic assembler that can be used for any standard DG scheme with flux and source.
Definition ddg.hh:192
SparsityPattern sparsity_pattern_jacobian
Definition ddg.hh:1076
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition ddg.hh:215
virtual void refinement_indicator(Vector< double > &indicator, const VectorType &solution_global) override
refinement indicator for adaptivity. Calls the model's cell_indicator and face_indicator functions.
Definition ddg.hh:266
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 ddg.hh:645
double average_time_jacobian_assembly()
Definition ddg.hh:1053
virtual const SparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition ddg.hh:254
std::vector< double > timings_jacobian
Definition ddg.hh:1081
typename Discretization::NumberType NumberType
Definition ddg.hh:198
std::vector< double > timings_residual
Definition ddg.hh:1080
virtual const SparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition ddg.hh:258
virtual void jacobian_mass(SparseMatrix< NumberType > &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha, NumberType beta) override
Definition ddg.hh:584
typename Discretization::VectorType VectorType
Definition ddg.hh:199
virtual void mass(VectorType &mass, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight) override
Definition ddg.hh:355
virtual void rebuild_jacobian_sparsity() override
Definition ddg.hh:242
std::vector< types::global_dof_index > extractor_dof_indices
Definition common.hh:355
typename Discretization::Components Components
Definition ddg.hh:201
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 ddg.hh:407
virtual void reinit_vector(VectorType &vec) const override
Definition ddg.hh:213
SparsityPattern sparsity_pattern_mass
Definition ddg.hh:1075
double average_time_residual_assembly()
Definition ddg.hh:1043
Assembler(Discretization &discretization, Model &model, const JSONValue &json)
Definition ddg.hh:204
Definition complex_math.hh:10
Definition tuples.hh:34
Definition ddg.hh:176
std::array< uint, 2 > cell_indices
Definition ddg.hh:177
std::array< double, 2 > values
Definition ddg.hh:178
Definition ddg.hh:175
std::vector< CopyFaceData_I > face_data
Definition ddg.hh:180
Definition ddg.hh:145
FullMatrix< NumberType > extractor_cell_jacobian
Definition ddg.hh:147
std::vector< types::global_dof_index > joint_dof_indices
Definition ddg.hh:148
FullMatrix< NumberType > cell_jacobian
Definition ddg.hh:146
void reinit(const FEInterfaceValues< dim > &fe_iv, uint n_extractors)
Definition ddg.hh:150
Definition ddg.hh:144
void reinit(const Iterator &cell, uint dofs_per_cell, uint n_extractors)
Definition ddg.hh:164
FullMatrix< NumberType > extractor_cell_jacobian
Definition ddg.hh:160
std::vector< CopyDataFace_J > face_data
Definition ddg.hh:162
FullMatrix< NumberType > cell_jacobian
Definition ddg.hh:159
std::vector< types::global_dof_index > local_dof_indices
Definition ddg.hh:161
Definition ddg.hh:119
std::vector< types::global_dof_index > joint_dof_indices
Definition ddg.hh:121
void reinit(const FEInterfaceValues< dim > &fe_iv)
Definition ddg.hh:123
Vector< NumberType > cell_residual
Definition ddg.hh:120
Definition ddg.hh:118
Vector< NumberType > cell_residual
Definition ddg.hh:130
std::vector< CopyDataFace_R > face_data
Definition ddg.hh:132
void reinit(const Iterator &cell, uint dofs_per_cell)
Definition ddg.hh:134
std::vector< types::global_dof_index > local_dof_indices
Definition ddg.hh:131
Class to hold data for each assembly thread, i.e. FEValues for cells, interfaces, as well as pre-allo...
Definition ddg.hh:37
std::vector< Vector< NumberType > > solution_dot
Definition ddg.hh:107
std::vector< Tensor< 1, dim > > cached_shape_grads
Definition ddg.hh:114
std::vector< double > cached_shape_values
Definition ddg.hh:113
std::vector< Vector< NumberType > > solution
Definition ddg.hh:104
typename Discretization::NumberType NumberType
Definition ddg.hh:39
std::vector< std::vector< Tensor< 1, dim, NumberType > > > solution_grad
Definition ddg.hh:105
array< std::vector< std::vector< Tensor< 1, dim, NumberType > > >, 2 > solution_grad_interface
Definition ddg.hh:109
ScratchData(const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< dim > &quadrature, const 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 ddg.hh:41
ScratchData(const ScratchData< Discretization > &scratch_data)
Definition ddg.hh:71
array< std::vector< Vector< NumberType > >, 2 > solution_interface
Definition ddg.hh:108
std::vector< Tensor< 2, dim > > cached_shape_hessians
Definition ddg.hh:115
FEInterfaceValues< dim > fe_interface_values
Definition ddg.hh:102
std::vector< std::vector< Tensor< 2, dim, NumberType > > > solution_hess
Definition ddg.hh:106
array< std::vector< std::vector< Tensor< 2, dim, NumberType > > >, 2 > solution_hess_interface
Definition ddg.hh:110
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