/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/discretization/FEM/assembler/ldg.hh Source File#
|
DiFfRG
|
ldg.hh
Go to the documentation of this file.
54 threads(json.get_uint("/discretization/threads")), batch_size(json.get_uint("/discretization/batch_size")),
101 virtual void refinement_indicator(Vector<double> & /*indicator*/, const VectorType & /*solution*/) = 0;
204 old_fe->get_mapping(), old_fe->get_fe(), old_fe->get_quadrature(), old_fe->get_update_flags()));
206 old_fe_i->get_mapping(), old_fe_i->get_fe(), old_fe_i->get_quadrature(), old_fe_i->get_update_flags()));
232 const auto &new_fe_interface_values(const t_Iterator &t_cell, uint f, uint sf, const t_Iterator &t_ncell,
320 CopyDataFace_J &new_face_data(const FEInterfaceValues<dim> &fe_iv_from, const FEInterfaceValues<dim> &fe_iv_to)
363 CopyDataFace_J &new_face_data(const array<unique_ptr<FEInterfaceValues<dim>>, n_fe_subsystems> &fe_iv,
373 copy_data_face.extractor_cell_jacobian.reinit(fe_iv[0]->n_current_interface_dofs(), n_extractors);
414 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "LDG1", "extractors", "variables">>(t);
416 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "LDG1", "LDG2", "extractors", "variables">>(
428 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "LDG1", "fe_derivatives", "fe_hessians",
431 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "LDG1", "LDG2", "fe_derivatives",
434 return named_tuple<std::tuple<T &...>, StringSet<"fe_functions", "LDG1", "LDG2", "LDG3", "fe_derivatives",
463 virtual void reinit_vector(VectorType &vec) const override { vec.reinit(dof_handler.n_dofs()); }
471 virtual void attach_data_output(DataOutput<dim, VectorType> &data_out, const VectorType &solution,
513 DoFTools::make_sparsity_pattern(*(dof_handler_list[0]), dsp, discretization.get_constraints(0), true);
534 build_ldg_sparsity(sparsity_pattern_jacobian, *(dof_handler_list[0]), *(dof_handler_list[0]), stencil,
545 // build the subjacobian sparsity patterns of all matrices that contribute to the jacobian = uu + ug*gu
550 build_ldg_sparsity(sparsity_pattern_wg[i], *(dof_handler_list[i]), *(dof_handler_list[i - 1]), 1);
574 build_ldg_sparsity(sparsity_pattern_jacobian, *(dof_handler_list[0]), *(dof_handler_list[0]), stencil, true);
579 virtual void refinement_indicator(Vector<double> &indicator, const VectorType &solution_global) override
585 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
608 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
637 model.face_indicator(local_indicator, normals[q_index], x_q, ref_conv(sol_q_s), ref_conv(sol_q_n));
657 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
665 virtual const BlockSparseMatrix<NumberType> &get_mass_matrix() const override { return mass_matrix; }
674 virtual void mass(VectorType &residual, const VectorType &solution_global, const VectorType &solution_global_dot,
682 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
719 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
730 virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight,
745 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
767 auto sol_q = std::tuple_cat(local_sol_q(solution, q_index), std::tie(extracted_data, variables));
785 const auto boundary_worker = [&](const Iterator &t_cell, const uint &face_no, Scratch &scratch_data,
804 auto sol_q = std::tuple_cat(local_sol_q(solution, q_index), std::tie(extracted_data, variables));
812 scalar_product(numflux[component_i], normals[q_index])); // phi_i(x_q) * numflux(x_q, u_q) * n(x_q)
816 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
848 auto sol_q_s = std::tuple_cat(local_sol_q(solution[0], q_index), std::tie(extracted_data, variables));
849 auto sol_q_n = std::tuple_cat(local_sol_q(solution[1], q_index), std::tie(extracted_data, variables));
880 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
886 virtual void jacobian_mass(BlockSparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
895 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
926 JxW[q_index] * fe_v[0]->shape_value_component(j, q_index, component_j) * // weight * dx * phi_j(x_q)
943 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
954 virtual void jacobian(BlockSparseMatrix<NumberType> &jacobian, const VectorType &solution_global,
976 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
997 SimpleMatrix<Tensor<1, dim, NumberType>, Components::count_fe_functions(), Components::count_extractors()>
999 SimpleMatrix<NumberType, Components::count_fe_functions(), Components::count_extractors()> j_extr_source;
1006 auto sol_q = std::tuple_cat(local_sol_q(solution, q_index), std::tie(extracted_data, variables));
1009 this->model.template jacobian_mass<0>(j_mass, x_q, solution[0][q_index], solution_dot[q_index]);
1010 this->model.template jacobian_mass<1>(j_mass_dot, x_q, solution[0][q_index], solution_dot[q_index]);
1059 const auto boundary_worker = [&](const Iterator &t_cell, const uint &face_no, Scratch &scratch_data,
1079 auto sol_q = std::tuple_cat(local_sol_q(solution, q_index), std::tie(extracted_data, variables));
1083 model.template jacobian_boundary_numflux<k, 0>(std::get<k>(j_boundary_numflux), normals[q_index], x_q,
1106 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
1152 auto sol_q_s = std::tuple_cat(local_sol_q(solution[0], q_index), std::tie(extracted_data, variables));
1153 auto sol_q_n = std::tuple_cat(local_sol_q(solution[1], q_index), std::tie(extracted_data, variables));
1155 model.template jacobian_numflux<k, 0>(std::get<k>(j_numflux), normals[q_index], x_q, fe_conv(sol_q_s),
1158 if (!std::get<k>(j_numflux)[0].is_finite() || !std::get<k>(j_numflux)[1].is_finite()) exception = true;
1181 constraints.distribute_local_to_global(c.cell_mass_jacobian, c.local_dof_indices[0], jacobian);
1183 constraints.distribute_local_to_global(cdf.cell_jacobian[0], cdf.joint_dof_indices[0], jacobian);
1193 FullMatrix<NumberType> extractor_dependence(c.local_dof_indices[0].size(), extractor_dof_indices.size());
1233 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
1239 tbb::blocked_range<uint>(1, Components::count_fe_subsystems()), [&](tbb::blocked_range<uint> rk) {
1244 tbb::blocked_range<uint>(0, Components::count_fe_functions(0)), [&](tbb::blocked_range<uint> r) {
1275 ss << " Reinit: " << average_time_reinit() * 1000 << "ms (" << num_reinits() << ")" << std::endl;
1276 ss << " Residual: " << average_time_residual_assembly() * 1000 << "ms (" << num_residuals() << ")"
1278 ss << " Jacobian: " << average_time_jacobian_assembly() * 1000 << "ms (" << num_jacobians() << ")"
1410 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
1439 const auto boundary_worker = [&](const Iterator &t_cell, const uint &face_no, Scratch &scratch_data,
1455 model.template ldg_boundary_numflux<to>(numflux, normals[q_index], x_q, solution[from][q_index]);
1462 scalar_product(numflux[component_i], normals[q_index])); // phi_i(x_q) * numflux(x_q, u_q) * n(x_q)
1466 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
1505 constraints.distribute_local_to_global(cdf.cell_residual, cdf.joint_dof_indices, ldg_vector_tmp);
1515 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
1532 void build_ldg_jacobian(const VectorType &solution_global, BlockSparseMatrix<NumberType> &ldg_jacobian,
1540 const auto cell_worker = [&](const Iterator &t_cell, Scratch &scratch_data, CopyData ©_data) {
1553 SimpleMatrix<Tensor<1, dim>, Components::count_fe_functions(to), Components::count_fe_functions(from)> j_flux;
1554 SimpleMatrix<NumberType, Components::count_fe_functions(to), Components::count_fe_functions(from)> j_source;
1574 const auto boundary_worker = [&](const Iterator &t_cell, const uint &face_no, Scratch &scratch_data,
1588 SimpleMatrix<Tensor<1, dim>, Components::count_fe_functions(to), Components::count_fe_functions(from)>
1609 const auto face_worker = [&](const Iterator &t_cell, const uint &f, const uint &sf, const Iterator &t_ncell,
1626 array<SimpleMatrix<Tensor<1, dim>, Components::count_fe_functions(to), Components::count_fe_functions(from)>,
1631 model.template jacobian_numflux<from, to>(j_numflux, normals[q_index], x_q, solution[0][from][q_index],
1637 const auto &component_i = fe_iv[to]->get_fe().system_to_component_index(cd_i[face_no_i]).first;
1641 const auto &component_j = fe_iv[from]->get_fe().system_to_component_index(cd_j[face_no_j]).first;
1668 MeshWorker::mesh_loop(dof_handler.begin_active(), dof_handler.end(), cell_worker, copier, scratch_data,
1671 component_mass_matrix_inverse.mmult(ldg_jacobian.block(c[0], c[1]), ldg_jacobian_tmp.block(c[0], c[1]),
1683 void build_ldg_sparsity(BlockSparsityPattern &sparsity_pattern, const DoFHandler<dim> &to_dofh,
1800 EoM_cell, solution_global, dof_handler, mapping, EoMfun, [&](const auto &p, const auto &) { return p; },
1810 update_values | update_gradients | update_quadrature_points | update_JxW_values | update_hessians));
1841 auto solution_tuple = std::tuple_cat(vector_to_tuple<Components::count_fe_subsystems()>(solutions_vector),
1849 void extract(std::array<NumberType, Components::count_extractors()> &data, const VectorType &solution_global,
1858 [&](const auto &p, const auto &values) { return postprocess ? model.EoM_postprocess(p, values) : p; },
1874 update_values | update_gradients | update_quadrature_points | update_JxW_values | update_hessians));
1900 auto solution_tuple = std::tuple_cat(vector_to_tuple<Components::count_fe_subsystems()>(solutions_vector),
1906 bool jacobian_extractors(FullMatrix<NumberType> &extractor_jacobian, const VectorType &solution_global,
1917 [&](const auto &p, const auto &values) { return model.EoM_postprocess(p, values); }, EoM_abs_tol,
1929 update_values | update_gradients | update_quadrature_points | update_JxW_values | update_hessians));
1963 auto solution_tuple = std::tuple_cat(vector_to_tuple<Components::count_fe_subsystems()>(solutions_vector),
1967 model.template jacobian_extractors<0>(extractor_jacobian_u, EoM, fe_more_conv(solution_tuple));
This is the general assembler interface for any kind of discretization. An assembler is responsible f...
Definition abstract_assembler.hh:39
Class to manage writing to files. FEM functions are written to vtk files and other data is written to...
Definition data_output.hh:20
FEOutput< dim, VectorType > & fe_output()
Returns a reference to the FEOutput object used to write FEM functions to .vtu files and ....
The LDG assembler that can be used for any LDG scheme, with as many levels as one wants.
Definition ldg.hh:397
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 ldg.hh:1849
double average_time_residual_assembly() const
Definition ldg.hh:1293
array< BlockVector< NumberType >, Components::count_fe_subsystems()> sol_vector
Definition ldg.hh:1327
void build_ldg_jacobian(const VectorType &solution_global, BlockSparseMatrix< NumberType > &ldg_jacobian, BlockSparseMatrix< NumberType > &ldg_jacobian_tmp) const
Build the LDG jacobian at level 'to', which takes information from level 'from'.
Definition ldg.hh:1532
SparseMatrix< NumberType > component_mass_matrix_inverse
Definition ldg.hh:1340
BlockSparsityPattern sparsity_pattern_jacobian
Definition ldg.hh:1331
virtual void residual(VectorType &residual, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType weight_mass, const VectorType &variables=VectorType()) override
Construct the system residual, i.e. Res = grad(flux) - source.
Definition ldg.hh:730
virtual void jacobian_mass(BlockSparseMatrix< NumberType > &jacobian, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType alpha=1., NumberType beta=1.) override
Definition ldg.hh:886
array< BlockSparseMatrix< NumberType >, Components::count_fe_subsystems()> j_wg_tmp
Definition ldg.hh:1344
virtual const BlockSparseMatrix< NumberType > & get_mass_matrix() const override
Obtain the mass matrix.
Definition ldg.hh:665
void build_ldg_sparsity(BlockSparsityPattern &sparsity_pattern, const DoFHandler< dim > &to_dofh, const DoFHandler< dim > &from_dofh, const int stencil=1, bool add_extractor_dofs=false) const
Create a sparsity pattern for matrices between the DoFs of two DoFHandlers, with given stencil size.
Definition ldg.hh:1683
double average_time_jacobian_assembly() const
Definition ldg.hh:1303
void readouts(DataOutput< dim, VectorType > &data_out, const VectorType &solution_global, const VectorType &variables) const
Definition ldg.hh:1794
typename Discretization::NumberType NumberType
Definition ldg.hh:403
std::vector< double > timings_residual
Definition ldg.hh:1347
array< BlockSparsityPattern, Components::count_fe_subsystems()> sparsity_pattern_ug
Definition ldg.hh:1333
array< BlockVector< NumberType >, Components::count_fe_subsystems()> sol_vector_tmp
Definition ldg.hh:1328
auto fe_more_conv(std::tuple< T &... > &t) const
Definition ldg.hh:425
virtual void residual_variables(VectorType &residual, const VectorType &variables, const VectorType &spatial_solution) override
Definition ldg.hh:1996
typename Discretization::VectorType VectorType
Definition ldg.hh:404
BlockSparsityPattern sparsity_pattern_mass
Definition ldg.hh:1332
virtual void reinit_vector(VectorType &vec) const override
Definition ldg.hh:463
virtual void jacobian_variables(FullMatrix< NumberType > &jacobian, const VectorType &variables, const VectorType &spatial_solution) override
Definition ldg.hh:2009
typename Discretization::Components Components
Definition ldg.hh:406
virtual void attach_data_output(DataOutput< dim, VectorType > &data_out, const VectorType &solution, const VectorType &variables, const VectorType &dt_solution=VectorType(), const VectorType &residual=VectorType()) override
Attach all intermediate (ldg) vectors to the data output.
Definition ldg.hh:471
virtual void mass(VectorType &residual, const VectorType &solution_global, const VectorType &solution_global_dot, NumberType weight) override
Construct the mass.
Definition ldg.hh:674
std::vector< const DoFHandler< dim > * > dof_handler_list
Definition ldg.hh:1325
array< BlockSparsityPattern, Components::count_fe_subsystems()> sparsity_pattern_wg
Definition ldg.hh:1335
void rebuild_ldg_jacobian(const VectorType &sol) const
Definition ldg.hh:1373
array< bool, Components::count_fe_subsystems()> ldg_matrix_built
Definition ldg.hh:1350
array< BlockSparseMatrix< NumberType >, Components::count_fe_subsystems()> j_ug
Definition ldg.hh:1341
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition ldg.hh:496
virtual void rebuild_jacobian_sparsity() override
Definition ldg.hh:572
bool jacobian_extractors(FullMatrix< NumberType > &extractor_jacobian, const VectorType &solution_global, const VectorType &variables)
Definition ldg.hh:1906
void build_inverse(const SparseMatrix< NumberType > &in, SparseMatrix< NumberType > &out) const
Build the inverse of matrix in and save the result to out.
Definition ldg.hh:1767
array< BlockSparsityPattern, Components::count_fe_subsystems()> sparsity_pattern_gu
Definition ldg.hh:1334
array< BlockSparseMatrix< NumberType >, Components::count_fe_subsystems()> jacobian_tmp
Definition ldg.hh:1337
virtual const BlockSparsityPattern & get_sparsity_pattern_jacobian() const override
Obtain the sparsity pattern of the jacobian matrix.
Definition ldg.hh:661
Assembler(Discretization &discretization, Model &model, const JSONValue &json)
Definition ldg.hh:453
BlockSparseMatrix< NumberType > mass_matrix
Definition ldg.hh:1339
virtual void jacobian(BlockSparseMatrix< NumberType > &jacobian, const VectorType &solution_global, NumberType weight, const VectorType &solution_global_dot, NumberType alpha, NumberType beta, const VectorType &variables=VectorType()) override
Construct the system jacobian, i.e. dRes/du.
Definition ldg.hh:954
array< BlockSparseMatrix< NumberType >, Components::count_fe_subsystems()> j_wg
Definition ldg.hh:1343
std::vector< types::global_dof_index > extractor_dof_indices
Definition ldg.hh:142
void rebuild_ldg_vectors(const VectorType &sol) const
Definition ldg.hh:1357
array< bool, Components::count_fe_subsystems()> jacobian_tmp_built
Definition ldg.hh:1351
array< Vector< NumberType >, Components::count_fe_subsystems()> sol_vector_vec_tmp
Definition ldg.hh:1329
std::vector< double > timings_jacobian
Definition ldg.hh:1348
array< BlockSparseMatrix< NumberType >, Components::count_fe_subsystems()> j_gu
Definition ldg.hh:1342
void build_ldg_vector(const VectorType &solution_global, VectorTypeldg &ldg_vector, VectorTypeldg &ldg_vector_tmp) const
Build the LDG vector at level 'to', which takes information from level 'from'.
Definition ldg.hh:1401
virtual void refinement_indicator(Vector< double > &indicator, const VectorType &solution_global) override
Definition ldg.hh:579
Definition ldg.hh:41
FullMatrix< NumberType > extractor_jacobian_du
Definition ldg.hh:140
uint num_variable_jacobians() const
Definition ldg.hh:121
FullMatrix< NumberType > extractor_jacobian_u
Definition ldg.hh:139
FullMatrix< NumberType > extractor_jacobian_ddu
Definition ldg.hh:141
typename Discretization::NumberType NumberType
Definition ldg.hh:45
double average_time_variable_residual_assembly()
Definition ldg.hh:103
FullMatrix< NumberType > extractor_jacobian
Definition ldg.hh:138
DoFHandler< dim >::cell_iterator old_EoM_cell
Definition ldg.hh:134
typename Discretization::VectorType VectorType
Definition ldg.hh:46
virtual void set_time(double t) override
Set the current time. The assembler should usually just forward this to the numerical model.
Definition ldg.hh:99
const auto & get_discretization() const
Definition ldg.hh:70
virtual void refinement_indicator(Vector< double > &, const VectorType &)=0
virtual void reinit() override
Reinitialize the assembler. This is necessary if the mesh has changed, e.g. after a mesh refinement.
Definition ldg.hh:73
double average_time_variable_jacobian_assembly()
Definition ldg.hh:113
virtual IndexSet get_differential_indices() const override
Obtain the dofs which contain time derivatives.
Definition ldg.hh:64
std::vector< double > timings_variable_residual
Definition ldg.hh:144
std::vector< double > timings_variable_jacobian
Definition ldg.hh:145
typename Discretization::Components Components
Definition ldg.hh:48
LDGAssemblerBase(Discretization &discretization, Model &model, const JSONValue &json)
Definition ldg.hh:51
DoFHandler< dim >::cell_iterator EoM_cell
Definition ldg.hh:133
std::vector< types::global_dof_index > extractor_dof_indices
Definition ldg.hh:142
uint num_variable_residuals() const
Definition ldg.hh:111
virtual void rebuild_jacobian_sparsity()=0
const DoFHandler< dim > & dof_handler
Definition ldg.hh:127
Definition quadrature.hh:56
size_t size() const
A simple NxM-matrix class, which is used for cell-wise Jacobians.
Definition tuples.hh:156
Definition complex_math.hh:10
dealii::Point< dim > get_EoM_point(typename dealii::DoFHandler< dim >::cell_iterator &EoM_cell, const VectorType &sol, const dealii::DoFHandler< dim > &dof_handler, const dealii::Mapping< dim > &mapping, const EoMFUN &get_EoM, const EoMPFUN &EoM_postprocess=[](const auto &p, const auto &values) { return p;}, const double EoM_abs_tol=1e-5, const uint max_iter=100)
Get the EoM point for a given solution and model.
Definition eom.hh:674
auto local_sol_q(const std::array< T, N > &a, uint q_index)
Definition tuples.hh:278
bool KOKKOS_INLINE_FUNCTION is_close(T1 a, T2 b, T3 eps_)
Function to evaluate whether two floats are equal to numerical precision. Tests for both relative and...
Definition math.hh:168
constexpr void constexpr_for(F &&f)
A compile-time for loop, which calls the lambda f of signature void(integer) for each index.
Definition utils.hh:31
Definition ldg.hh:379
std::array< uint, 2 > cell_indices
Definition ldg.hh:380
std::array< double, 2 > values
Definition ldg.hh:381
Definition ldg.hh:378
std::vector< CopyFaceData_I > face_data
Definition ldg.hh:383
Definition ldg.hh:298
FullMatrix< NumberType > cell_jacobian
Definition ldg.hh:299
std::vector< types::global_dof_index > joint_dof_indices_from
Definition ldg.hh:300
std::vector< types::global_dof_index > joint_dof_indices_to
Definition ldg.hh:301
Definition ldg.hh:332
FullMatrix< NumberType > extractor_cell_jacobian
Definition ldg.hh:334
array< FullMatrix< NumberType >, n_fe_subsystems > cell_jacobian
Definition ldg.hh:333
array< vector< types::global_dof_index >, n_fe_subsystems > joint_dof_indices
Definition ldg.hh:335
Definition ldg.hh:331
void reinit(const array< Iterator, n_fe_subsystems > &cell, const uint n_extractors)
Definition ldg.hh:346
CopyDataFace_J & new_face_data(const array< unique_ptr< FEInterfaceValues< dim > >, n_fe_subsystems > &fe_iv, const uint n_extractors)
Definition ldg.hh:363
FullMatrix< NumberType > cell_mass_jacobian
Definition ldg.hh:339
array< FullMatrix< NumberType >, n_fe_subsystems > cell_jacobian
Definition ldg.hh:338
vector< CopyDataFace_J > face_data
Definition ldg.hh:342
FullMatrix< NumberType > extractor_cell_jacobian
Definition ldg.hh:340
array< vector< types::global_dof_index >, n_fe_subsystems > local_dof_indices
Definition ldg.hh:341
Definition ldg.hh:297
CopyDataFace_J & new_face_data(const FEInterfaceValues< dim > &fe_iv_from, const FEInterfaceValues< dim > &fe_iv_to)
Definition ldg.hh:320
std::vector< types::global_dof_index > local_dof_indices_to
Definition ldg.hh:306
FullMatrix< NumberType > cell_jacobian
Definition ldg.hh:304
std::vector< CopyDataFace_J > face_data
Definition ldg.hh:307
std::vector< types::global_dof_index > local_dof_indices_from
Definition ldg.hh:305
void reinit(const Iterator &cell_from, const Iterator &cell_to, uint dofs_per_cell_from, uint dofs_per_cell_to)
Definition ldg.hh:310
Definition ldg.hh:267
Vector< NumberType > cell_residual
Definition ldg.hh:268
std::vector< types::global_dof_index > joint_dof_indices
Definition ldg.hh:269
Definition ldg.hh:266
std::vector< CopyDataFace_R > face_data
Definition ldg.hh:275
void reinit(const Iterator &cell, uint dofs_per_cell)
Definition ldg.hh:277
std::vector< types::global_dof_index > local_dof_indices
Definition ldg.hh:274
Vector< NumberType > cell_residual
Definition ldg.hh:272
CopyDataFace_R & new_face_data(const FEInterfaceValues< dim > &fe_iv)
Definition ldg.hh:287
Class to hold data for each assembly thread, i.e. FEValues for cells, interfaces, as well as pre-allo...
Definition ldg.hh:154
array< uint, n_fe_subsystems > n_components
Definition ldg.hh:251
typename Discretization::NumberType NumberType
Definition ldg.hh:156
vector< Vector< NumberType > > solution_dot
Definition ldg.hh:262
static constexpr uint n_fe_subsystems
Definition ldg.hh:157
const auto & new_fe_values(const t_Iterator &t_cell)
Definition ldg.hh:224
array< Iterator, n_fe_subsystems > cell
Definition ldg.hh:252
array< unique_ptr< FEFaceValues< dim > >, n_fe_subsystems > fe_boundary_values
Definition ldg.hh:257
typename Triangulation< dim >::active_cell_iterator t_Iterator
Definition ldg.hh:159
ScratchData(const ScratchData< Discretization > &scratch_data)
Definition ldg.hh:196
typename DoFHandler< dim >::active_cell_iterator Iterator
Definition ldg.hh:158
array< unique_ptr< FEInterfaceValues< dim > >, n_fe_subsystems > fe_interface_values
Definition ldg.hh:256
array< array< vector< Vector< NumberType > >, n_fe_subsystems >, 2 > solution_interface
Definition ldg.hh:263
ScratchData(const Mapping< dim > &mapping, const vector< const DoFHandler< dim > * > &dofh, const Quadrature< dim > &quadrature, const Quadrature< dim - 1 > &quadrature_face, const UpdateFlags update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values, const UpdateFlags interface_update_flags=update_values|update_gradients|update_quadrature_points|update_JxW_values|update_normal_vectors)
Definition ldg.hh:161
const auto & new_fe_interface_values(const t_Iterator &t_cell, uint f, uint sf, const t_Iterator &t_ncell, uint nf, unsigned int nsf)
Definition ldg.hh:232
array< Iterator, n_fe_subsystems > ncell
Definition ldg.hh:253
array< std::vector< uint >, n_fe_subsystems > comp
Definition ldg.hh:259
array< unique_ptr< FEValues< dim > >, n_fe_subsystems > fe_values
Definition ldg.hh:255
array< vector< Vector< NumberType > >, n_fe_subsystems > solution
Definition ldg.hh:261
const auto & new_fe_boundary_values(const t_Iterator &t_cell, uint face_no)
Definition ldg.hh:242
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