17 #include <deal.II/fe/fe_enriched.h> 18 #include <deal.II/fe/fe_tools.h> 20 #include <deal.II/base/std_cxx14/memory.h> 23 DEAL_II_NAMESPACE_OPEN
35 std::vector<unsigned int>
36 build_multiplicities(
const std::vector<std::vector<T > > &functions )
38 std::vector<unsigned int> multiplicities;
39 multiplicities.push_back(1);
40 for (
unsigned int i = 0; i < functions.size(); i++)
41 multiplicities.push_back(functions[i].size());
43 return multiplicities;
50 template <
int dim,
int spacedim>
51 std::vector< const FiniteElement< dim, spacedim > * >
55 std::vector< const FiniteElement< dim, spacedim > * > fes;
56 fes.push_back(fe_base);
57 for (
unsigned int i = 0; i < fe_enriched.size(); i++)
58 fes.push_back(fe_enriched[i]);
68 template <
int dim,
int spacedim>
71 const std::vector< unsigned int > &multiplicities,
72 const std::vector<std::vector<std::function<
const Function<spacedim> *(
const typename ::Triangulation<dim, spacedim>::cell_iterator &) > > > &functions)
77 ExcMessage(
"FEs and multiplicities should have the same size"));
85 const unsigned int n_comp_base = fes[0]->n_components();
88 for (
unsigned int fe=1; fe < fes.size(); fe++)
93 ExcMessage(
"Only dominating FE_Nothing can be used in FE_Enriched"));
96 ExcMessage(
"All elements must have the same number of components"));
105 template <
int dim,
int spacedim>
110 for (
unsigned int fe=1; fe < fes.size(); fe++)
122 template <
int dim,
int spacedim>
132 template <
int dim,
int spacedim>
146 return enrichment_function;
153 template <
int dim,
int spacedim>
164 template <
int dim,
int spacedim>
166 const std::vector< unsigned int > &multiplicities,
169 FiniteElement<dim,spacedim> (
FETools::Compositing::multiply_dof_numbers(fes,multiplicities,false),
170 FETools::Compositing::compute_restriction_is_additive_flags(fes,multiplicities),
171 FETools::Compositing::compute_nonzero_components(fes,multiplicities,false)),
172 enrichments(functions),
177 Assert(internal::FE_Enriched::consistency_check(fes,multiplicities,functions),
185 for (
unsigned int fe=1; fe < fes.size(); fe++)
193 for (
unsigned int system_index=0; system_index<this->
dofs_per_cell;
203 ExcMessage(
"Size mismatch for base_no_mult_local_enriched_dofs: " 205 "; base_no = " + std::to_string(base_no) +
206 "; base_m = " + std::to_string(base_m) +
207 "; system_index = " + std::to_string(system_index)));
228 template <
int dim,
int spacedim>
236 template <
int dim,
int spacedim>
241 Assert(!is_enriched,
ExcMessage(
"For enriched finite elements shape_value() can not be defined on the reference element."));
242 return fe_system->shape_value(i,p);
246 template <
int dim,
int spacedim>
247 std::unique_ptr<FiniteElement<dim,spacedim> >
250 std::vector< const FiniteElement< dim, spacedim > * > fes;
251 std::vector< unsigned int > multiplicities;
253 for (
unsigned int i=0; i<this->n_base_elements(); i++)
255 fes.push_back( & base_element(i) );
256 multiplicities.push_back(this->element_multiplicity(i) );
265 template <
int dim,
int spacedim>
269 UpdateFlags out = fe_system->requires_update_flags(flags);
289 template <
int dim,
int spacedim>
291 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
298 auto update_each_flags = fes_data->update_each;
299 auto data = std_cxx14::make_unique<InternalData>(std::move(fes_data));
302 data->update_each = update_each_flags;
305 data->enrichment.resize(this->n_base_elements());
307 const unsigned int n_q_points = quadrature.
size();
309 for (
unsigned int base=0; base < this->n_base_elements(); ++base)
311 data->enrichment[base].resize(this->element_multiplicity(base));
312 for (
unsigned int m = 0; m < this->element_multiplicity(base); ++m)
315 data->enrichment[base][m].values.resize(n_q_points);
318 data->enrichment[base][m].gradients.resize(n_q_points);
321 data->enrichment[base][m].hessians.resize(n_q_points);
325 return std::move(data);
329 template <
int dim,
int spacedim>
330 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
336 auto data = fe_system->get_face_data(update_flags,mapping,quadrature,output_data);
343 template <
int dim,
int spacedim>
344 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
350 auto data = fe_system->get_subface_data(update_flags,mapping,quadrature,output_data);
357 template <
int dim,
int spacedim>
358 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
364 auto data = fe_system->get_data(flags,mapping,quadrature,output_data);
371 template <
int dim,
int spacedim>
373 const std::vector<unsigned int> &multiplicities)
375 Assert (fes.size() == multiplicities.size(),
379 this->base_to_block_indices.reinit(0, 0);
381 for (
unsigned int i=0; i<fes.size(); i++)
382 if (multiplicities[i]>0)
383 this->base_to_block_indices.push_back( multiplicities[i] );
388 this->system_to_component_table.resize(this->dofs_per_cell);
389 this->face_system_to_component_table.resize(this->dofs_per_face);
392 this->system_to_component_table,
393 this->component_to_base_table,
398 this->face_system_to_component_table,
407 this->interface_constraints = fe_system->interface_constraints;
416 this->unit_support_points = fe_system->unit_support_points;
417 this->unit_face_support_points = fe_system->unit_face_support_points;
422 this->adjust_line_dof_index_for_line_orientation_table = fe_system->adjust_line_dof_index_for_line_orientation_table;
427 template <
int dim,
int spacedim>
431 std::ostringstream namebuf;
433 namebuf <<
"FE_Enriched<" 436 for (
unsigned int i=0; i< this->n_base_elements(); ++i)
438 namebuf << base_element(i).get_name();
439 if (this->element_multiplicity(i) != 1)
440 namebuf <<
'^' << this->element_multiplicity(i);
441 if (i != this->n_base_elements()-1)
446 return namebuf.str();
450 template <
int dim,
int spacedim>
454 return fe_system->base_element(index);
458 template <
int dim,
int spacedim>
465 const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data,
470 Assert (dynamic_cast<const InternalData *> (&fe_internal) !=
nullptr,
472 const InternalData &fe_data =
static_cast<const InternalData &
> (fe_internal);
475 fe_system->fill_fe_values(cell,
481 *fe_data.fesystem_data,
485 multiply_by_enrichment(quadrature,
493 template <
int dim,
int spacedim>
497 const unsigned int face_no,
501 const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data,
506 Assert (dynamic_cast<const InternalData *> (&fe_internal) !=
nullptr,
508 const InternalData &fe_data =
static_cast<const InternalData &
> (fe_internal);
511 fe_system->fill_fe_face_values(cell,
517 *fe_data.fesystem_data,
521 multiply_by_enrichment(quadrature,
529 template <
int dim,
int spacedim>
533 const unsigned int face_no,
534 const unsigned int sub_no,
538 const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data,
543 Assert (dynamic_cast<const InternalData *> (&fe_internal) !=
nullptr,
545 const InternalData &fe_data =
static_cast<const InternalData &
> (fe_internal);
548 fe_system->fill_fe_subface_values(cell,
555 *fe_data.fesystem_data,
559 multiply_by_enrichment(quadrature,
568 template <
int dim,
int spacedim>
588 const unsigned int n_q_points = quadrature.
size();
598 for (
unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
601 base_fe = base_element(base_no);
609 for (
unsigned int system_index=0; system_index<this->dofs_per_cell;
611 if (this->system_to_base_table[system_index].first.first == base_no)
614 base_index = this->system_to_base_table[system_index].second;
623 unsigned int out_index = 0;
624 for (
unsigned int i=0; i<system_index; ++i)
625 out_index += this->n_nonzero_components(i);
626 unsigned int in_index = 0;
627 for (
unsigned int i=0; i<base_index; ++i)
631 Assert (this->n_nonzero_components(system_index) ==
634 for (
unsigned int s=0; s<this->n_nonzero_components(system_index); ++s)
637 for (
unsigned int q=0; q<n_q_points; ++q)
642 for (
unsigned int q=0; q<n_q_points; ++q)
647 for (
unsigned int q=0; q<n_q_points; ++q)
654 Assert (base_no_mult_local_enriched_dofs.size() == fe_data.
enrichment.size(),
658 for (
unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
660 Assert (base_no_mult_local_enriched_dofs[base_no].size() == fe_data.
enrichment[base_no].size(),
663 for (
unsigned int m=0; m < base_no_mult_local_enriched_dofs[base_no].size(); m++)
667 if (base_no_mult_local_enriched_dofs[base_no][m].size()==0)
670 Assert (enrichments[base_no-1][m](cell) !=
nullptr,
671 ExcMessage(
"The pointer to the enrichment function is NULL"));
674 ExcMessage(
"Only scalar-valued enrichment functions are allowed"));
681 for (
unsigned int q=0; q<n_q_points; q++)
690 for (
unsigned int q=0; q<n_q_points; q++)
699 for (
unsigned int q=0; q<n_q_points; q++)
712 for (
unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
714 for (
unsigned int m=0; m < base_no_mult_local_enriched_dofs[base_no].size(); m++)
715 for (
unsigned int i=0; i < base_no_mult_local_enriched_dofs[base_no][m].size(); i++)
717 const unsigned int enriched_dof = base_no_mult_local_enriched_dofs[base_no][m][i];
718 for (
unsigned int q=0; q<n_q_points; ++q)
734 for (
unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
736 for (
unsigned int m=0; m < base_no_mult_local_enriched_dofs[base_no].size(); m++)
737 for (
unsigned int i=0; i < base_no_mult_local_enriched_dofs[base_no][m].size(); i++)
739 const unsigned int enriched_dof = base_no_mult_local_enriched_dofs[base_no][m][i];
740 for (
unsigned int q=0; q<n_q_points; ++q)
750 for (
unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
752 for (
unsigned int m=0; m < base_no_mult_local_enriched_dofs[base_no].size(); m++)
753 for (
unsigned int i=0; i < base_no_mult_local_enriched_dofs[base_no][m].size(); i++)
755 const unsigned int enriched_dof = base_no_mult_local_enriched_dofs[base_no][m][i];
756 for (
unsigned int q=0; q<n_q_points; ++q)
765 template <
int dim,
int spacedim>
774 template <
int dim,
int spacedim>
783 template <
int dim,
int spacedim>
792 fe_system->get_face_interpolation_matrix(fe_enr_other->get_fe_system(),
798 ExcInterpolationNotImplemented()));
803 template <
int dim,
int spacedim>
807 const unsigned int subface,
813 fe_system->get_subface_interpolation_matrix(fe_enr_other->get_fe_system(),
820 ExcInterpolationNotImplemented()));
825 template <
int dim,
int spacedim>
826 std::vector<std::pair<unsigned int, unsigned int> >
833 return fe_system->hp_vertex_dof_identities(fe_enr_other->get_fe_system());
838 return std::vector<std::pair<unsigned int, unsigned int> >();
843 template <
int dim,
int spacedim>
844 std::vector<std::pair<unsigned int, unsigned int> >
851 return fe_system->hp_line_dof_identities(fe_enr_other->get_fe_system());
856 return std::vector<std::pair<unsigned int, unsigned int> >();
861 template <
int dim,
int spacedim>
862 std::vector<std::pair<unsigned int, unsigned int> >
869 return fe_system->hp_quad_dof_identities(fe_enr_other->get_fe_system());
874 return std::vector<std::pair<unsigned int, unsigned int> >();
879 template <
int dim,
int spacedim>
897 return fe_system->compare_for_face_domination(fe_enr_other->get_fe_system());
906 template <
int dim,
int spacedim>
911 return fe_system->get_prolongation_matrix(child, refinement_case);
915 template <
int dim,
int spacedim>
920 return fe_system->get_restriction_matrix(child, refinement_case);
927 template <
int dim,
int spacedim>
930 fesystem_data(
std::move(fesystem_data))
934 template <
int dim,
int spacedim>
939 return fesystem_data->get_fe_data(base_no);
943 template <
int dim,
int spacedim>
948 return fesystem_data->get_fe_output_object(base_no);
953 #include "fe_enriched.inst" 955 DEAL_II_NAMESPACE_CLOSE
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
bool is_dominating() const
InternalData(std::unique_ptr< typename FESystem< dim, spacedim >::InternalData > fesystem_data)
unsigned int n_nonzero_components(const unsigned int i) const
FE_Enriched(const FiniteElement< dim, spacedim > &fe_base, const FiniteElement< dim, spacedim > &fe_enriched, const Function< spacedim > *enrichment_function)
void initialize(const std::vector< const FiniteElement< dim, spacedim > *> &fes, const std::vector< unsigned int > &multiplicities)
SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > & get_fe_output_object(const unsigned int base_no) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
std::vector< std::vector< EnrichmentValues > > enrichment
Transformed quadrature points.
#define AssertThrow(cond, exc)
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
void multiply_by_enrichment(const Quadrature< dim_1 > &quadrature, const InternalData &fe_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename Triangulation< dim, spacedim >::cell_iterator &cell, internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
std::unique_ptr< To > dynamic_unique_cast(std::unique_ptr< From > &&p)
const FESystem< dim, spacedim > & get_fe_system() const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual bool hp_constraints_are_implemented() const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Third derivatives of shape functions.
#define Assert(cond, exc)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
FiniteElement< dim, spacedim >::InternalDataBase & get_fe_data(const unsigned int base_no) const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
virtual std::string get_name() const
SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
Second derivatives of shape functions.
std::string dim_string(const int dim, const int spacedim)
unsigned int size() const
const unsigned int dofs_per_cell
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Shape function gradients.
std::vector< std::vector< std::vector< unsigned int > > > base_no_mult_local_enriched_dofs
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
static ::ExceptionBase & ExcNotImplemented()
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > setup_data(std::unique_ptr< typename FESystem< dim, spacedim >::InternalData > fes_data, const UpdateFlags flags, const Quadrature< dim_1 > &quadrature) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
unsigned int n_base_elements() const
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
static ::ExceptionBase & ExcInternalError()
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const