Reference documentation for deal.II version 9.0.0
fe_enriched.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/fe/fe_enriched.h>
18 #include <deal.II/fe/fe_tools.h>
19 
20 #include <deal.II/base/std_cxx14/memory.h>
21 
22 
23 DEAL_II_NAMESPACE_OPEN
24 
25 namespace internal
26 {
27  namespace FE_Enriched
28  {
29  namespace
30  {
34  template <typename T>
35  std::vector<unsigned int>
36  build_multiplicities(const std::vector<std::vector<T > > &functions )
37  {
38  std::vector<unsigned int> multiplicities;
39  multiplicities.push_back(1); // the first one is non-enriched FE
40  for (unsigned int i = 0; i < functions.size(); i++)
41  multiplicities.push_back(functions[i].size());
42 
43  return multiplicities;
44  }
45 
46 
50  template <int dim, int spacedim>
51  std::vector< const FiniteElement< dim, spacedim > * >
52  build_fes(const FiniteElement<dim,spacedim> *fe_base,
53  const std::vector<const FiniteElement<dim,spacedim> * > &fe_enriched)
54  {
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]);
59 
60  return fes;
61  }
62 
63 
68  template <int dim, int spacedim>
69  bool
70  consistency_check (const std::vector< const FiniteElement< dim, spacedim > * > &fes,
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)
73  {
74  AssertThrow(fes.size() > 0,
75  ExcMessage("FEs size should be >=1"));
76  AssertThrow(fes.size() == multiplicities.size(),
77  ExcMessage("FEs and multiplicities should have the same size"));
78 
79  AssertThrow (functions.size() == fes.size() - 1,
80  ExcDimensionMismatch(functions.size(), fes.size()-1));
81 
82  AssertThrow(multiplicities[0] == 1,
83  ExcMessage("First multiplicity should be 1"));
84 
85  const unsigned int n_comp_base = fes[0]->n_components();
86 
87  // start from fe=1 as 0th is always non-enriched FE.
88  for (unsigned int fe=1; fe < fes.size(); fe++)
89  {
90  const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(fes[fe]);
91  if (fe_nothing)
92  AssertThrow (fe_nothing->is_dominating(),
93  ExcMessage("Only dominating FE_Nothing can be used in FE_Enriched"));
94 
95  AssertThrow (fes[fe]->n_components() == n_comp_base,
96  ExcMessage("All elements must have the same number of components"));
97  }
98  return true;
99  }
100 
101 
105  template <int dim, int spacedim>
106  bool
107  check_if_enriched (const std::vector< const FiniteElement< dim, spacedim > * > &fes)
108  {
109  // start from fe=1 as 0th is always non-enriched FE.
110  for (unsigned int fe=1; fe < fes.size(); fe++)
111  if (dynamic_cast<const FE_Nothing<dim>*>(fes[fe]) == nullptr)
112  // this is not FE_Nothing => there will be enrichment
113  return true;
114 
115  return false;
116  }
117  }
118  }
119 }
120 
121 
122 template <int dim, int spacedim>
124  :
125  FE_Enriched<dim,spacedim>(fe_base,
126  FE_Nothing<dim,spacedim>(fe_base.n_components(),true),
127  nullptr)
128 {
129 }
130 
131 
132 template <int dim, int spacedim>
134  const FiniteElement<dim,spacedim> &fe_enriched,
135  const Function<spacedim> *enrichment_function)
136  :
137  FE_Enriched<dim,spacedim>
138  (&fe_base,
139  std::vector<const FiniteElement<dim,spacedim>*>(1, &fe_enriched),
140  std::vector<std::vector<std::function<const Function<spacedim> *(const typename Triangulation<dim, spacedim>::cell_iterator &) > > >
141  (1,
142  std::vector<std::function<const Function<spacedim> *(const typename Triangulation<dim, spacedim>::cell_iterator &) > >
143  (1,
144  [=] (const typename Triangulation<dim, spacedim>::cell_iterator &) -> const Function<spacedim> *
145 {
146  return enrichment_function;
147 })
148  )
149 )
150 {}
151 
152 
153 template <int dim, int spacedim>
155  const std::vector<const FiniteElement<dim,spacedim> * > &fe_enriched,
156  const std::vector<std::vector<std::function<const Function<spacedim> *(const typename Triangulation<dim, spacedim>::cell_iterator &) > > > &functions)
157  :
158  FE_Enriched<dim,spacedim> (internal::FE_Enriched::build_fes(fe_base,fe_enriched),
159  internal::FE_Enriched::build_multiplicities(functions),
160  functions)
161 {}
162 
163 
164 template <int dim, int spacedim>
166  const std::vector< unsigned int > &multiplicities,
167  const std::vector<std::vector<std::function<const Function<spacedim> *(const typename Triangulation<dim, spacedim>::cell_iterator &) > > > &functions)
168  :
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),
173  is_enriched(internal::FE_Enriched::check_if_enriched(fes)),
174  fe_system(std_cxx14::make_unique<FESystem<dim, spacedim> >(fes,multiplicities))
175 {
176  // descriptive error are thrown within the function.
177  Assert(internal::FE_Enriched::consistency_check(fes,multiplicities,functions),
178  ExcInternalError());
179 
180  initialize(fes, multiplicities);
181 
182  // resize to be consistent with all FEs used to construct the FE_Enriched,
183  // even though we will never use the 0th element.
184  base_no_mult_local_enriched_dofs.resize(fes.size());
185  for (unsigned int fe=1; fe < fes.size(); fe++)
186  base_no_mult_local_enriched_dofs[fe].resize(multiplicities[fe]);
187 
190  this->n_base_elements()));
191 
192  // build the map: (base_no, base_m) -> vector of local element DoFs
193  for (unsigned int system_index=0; system_index<this->dofs_per_cell;
194  ++system_index)
195  {
196  const unsigned int base_no = this->system_to_base_table[system_index].first.first;
197  if (base_no == 0) // 0th is always non-enriched FE
198  continue;
199 
200  const unsigned int base_m = this->system_to_base_table[system_index].first.second;
201 
202  Assert (base_m < base_no_mult_local_enriched_dofs[base_no].size(),
203  ExcMessage("Size mismatch for base_no_mult_local_enriched_dofs: "
204  "base_index = " + std::to_string(this->system_to_base_table[system_index].second) +
205  "; base_no = " + std::to_string(base_no) +
206  "; base_m = " + std::to_string(base_m) +
207  "; system_index = " + std::to_string(system_index)));
208 
209  Assert (base_m < base_no_mult_local_enriched_dofs[base_no].size(),
210  ExcDimensionMismatch(base_m,
211  base_no_mult_local_enriched_dofs[base_no].size()));
212 
213  base_no_mult_local_enriched_dofs[base_no][base_m].push_back(system_index);
214  }
215 
216  // make sure that local_enriched_dofs.size() is correct, that is equals to DoFs
217  // per cell of the corresponding FE.
218  for (unsigned int base_no = 1; base_no < base_no_mult_local_enriched_dofs.size(); base_no++)
219  {
220  for (unsigned int m=0; m < base_no_mult_local_enriched_dofs[base_no].size(); m++)
221  Assert ( base_no_mult_local_enriched_dofs[base_no][m].size() == fes[base_no]->dofs_per_cell,
223  fes[base_no]->dofs_per_cell));
224  }
225 }
226 
227 
228 template <int dim, int spacedim>
229 const std::vector<std::vector<std::function<const Function<spacedim> *(const typename Triangulation<dim, spacedim>::cell_iterator &) > > >
231 {
232  return enrichments;
233 }
234 
235 
236 template <int dim, int spacedim>
237 double
239  const Point< dim > &p) const
240 {
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);
243 }
244 
245 
246 template <int dim, int spacedim>
247 std::unique_ptr<FiniteElement<dim,spacedim> >
249 {
250  std::vector< const FiniteElement< dim, spacedim > * > fes;
251  std::vector< unsigned int > multiplicities;
252 
253  for (unsigned int i=0; i<this->n_base_elements(); i++)
254  {
255  fes.push_back( & base_element(i) );
256  multiplicities.push_back(this->element_multiplicity(i) );
257  }
258 
259  return std::unique_ptr<FE_Enriched<dim,spacedim>>(new FE_Enriched<dim,spacedim>(fes,
260  multiplicities,
261  get_enrichments()));
262 }
263 
264 
265 template <int dim, int spacedim>
268 {
269  UpdateFlags out = fe_system->requires_update_flags(flags);
270 
271  if (is_enriched)
272  {
273  // if we ask for values or gradients, then we would need quadrature points
274  if (flags & (update_values | update_gradients))
276 
277  // if need gradients, add update_values due to product rule
278  if (out & update_gradients)
279  out |= update_values;
280  }
281 
282  Assert (!(flags & update_3rd_derivatives),
284 
285  return out;
286 }
287 
288 
289 template <int dim, int spacedim>
290 template <int dim_1>
291 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
293  const UpdateFlags flags,
294  const Quadrature<dim_1> &quadrature) const
295 {
296  // Pass ownership of the FiniteElement::InternalDataBase object
297  // that fes_data points to, to the new InternalData object.
298  auto update_each_flags = fes_data->update_each;
299  auto data = std_cxx14::make_unique<InternalData>(std::move(fes_data));
300 
301  // copy update_each from FESystem data:
302  data->update_each = update_each_flags;
303 
304  // resize cache array according to requested flags
305  data->enrichment.resize(this->n_base_elements());
306 
307  const unsigned int n_q_points = quadrature.size();
308 
309  for (unsigned int base=0; base < this->n_base_elements(); ++base)
310  {
311  data->enrichment[base].resize(this->element_multiplicity(base));
312  for (unsigned int m = 0; m < this->element_multiplicity(base); ++m)
313  {
314  if (flags & update_values)
315  data->enrichment[base][m].values.resize(n_q_points);
316 
317  if (flags & update_gradients)
318  data->enrichment[base][m].gradients.resize(n_q_points);
319 
320  if (flags & update_hessians)
321  data->enrichment[base][m].hessians.resize(n_q_points);
322  }
323  }
324 
325  return std::move(data);
326 }
327 
328 
329 template <int dim, int spacedim>
330 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
332  const Mapping<dim,spacedim> &mapping,
333  const Quadrature<dim-1> &quadrature,
335 {
336  auto data = fe_system->get_face_data(update_flags,mapping,quadrature,output_data);
337  return setup_data(Utilities::dynamic_unique_cast<typename FESystem<dim,spacedim>::InternalData>(std::move(data)),
338  update_flags,
339  quadrature);
340 }
341 
342 
343 template <int dim, int spacedim>
344 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
346  const Mapping<dim,spacedim> &mapping,
347  const Quadrature<dim-1> &quadrature,
349 {
350  auto data = fe_system->get_subface_data(update_flags,mapping,quadrature,output_data);
351  return setup_data(Utilities::dynamic_unique_cast<typename FESystem<dim,spacedim>::InternalData>(std::move(data)),
352  update_flags,
353  quadrature);
354 }
355 
356 
357 template <int dim, int spacedim>
358 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
360  const Mapping<dim,spacedim> &mapping,
361  const Quadrature<dim> &quadrature,
363 {
364  auto data = fe_system->get_data(flags,mapping,quadrature,output_data);
365  return setup_data(Utilities::dynamic_unique_cast<typename FESystem<dim,spacedim>::InternalData>(std::move(data)),
366  flags,
367  quadrature);
368 }
369 
370 
371 template <int dim, int spacedim>
373  const std::vector<unsigned int> &multiplicities)
374 {
375  Assert (fes.size() == multiplicities.size(),
376  ExcDimensionMismatch (fes.size(), multiplicities.size()) );
377 
378  // Note that we need to skip every fe with multiplicity 0 in the following block of code
379  this->base_to_block_indices.reinit(0, 0);
380 
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] );
384 
385  {
386  // If the system is not primitive, these have not been initialized by
387  // FiniteElement
388  this->system_to_component_table.resize(this->dofs_per_cell);
389  this->face_system_to_component_table.resize(this->dofs_per_face);
390 
391  FETools::Compositing::build_cell_tables(this->system_to_base_table,
392  this->system_to_component_table,
393  this->component_to_base_table,
394  *this,
395  false);
396 
397  FETools::Compositing::build_face_tables(this->face_system_to_base_table,
398  this->face_system_to_component_table,
399  *this,
400  false);
401  }
402 
403  // restriction and prolongation matrices are built on demand
404 
405  // now set up the interface constraints for h-refinement.
406  // take them from fe_system:
407  this->interface_constraints = fe_system->interface_constraints;
408 
409  // if we just wrap another FE (i.e. use FE_Nothing as a second FE)
410  // then it makes sense to have support points.
411  // However, functions like interpolate_boundary_values() need all FEs inside
412  // FECollection to be able to provide support points irrespectively whether
413  // this FE sits on the boundary or not. Thus for moment just copy support
414  // points from fe system:
415  {
416  this->unit_support_points = fe_system->unit_support_points;
417  this->unit_face_support_points = fe_system->unit_face_support_points;
418  }
419 
420  // take adjust_quad_dof_index_for_face_orientation_table from FESystem:
421  {
422  this->adjust_line_dof_index_for_line_orientation_table = fe_system->adjust_line_dof_index_for_line_orientation_table;
423  }
424 }
425 
426 
427 template <int dim, int spacedim>
428 std::string
430 {
431  std::ostringstream namebuf;
432 
433  namebuf << "FE_Enriched<"
434  << Utilities::dim_string(dim,spacedim)
435  << ">[";
436  for (unsigned int i=0; i< this->n_base_elements(); ++i)
437  {
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)
442  namebuf << '-';
443  }
444  namebuf << ']';
445 
446  return namebuf.str();
447 }
448 
449 
450 template <int dim, int spacedim>
452 FE_Enriched<dim,spacedim>::base_element (const unsigned int index) const
453 {
454  return fe_system->base_element(index);
455 }
456 
457 
458 template <int dim, int spacedim>
459 void
461  const CellSimilarity::Similarity cell_similarity,
462  const Quadrature< dim > &quadrature,
463  const Mapping< dim, spacedim > &mapping,
464  const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal,
465  const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data,
466  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
468  ) const
469 {
470  Assert (dynamic_cast<const InternalData *> (&fe_internal) != nullptr,
471  ExcInternalError());
472  const InternalData &fe_data = static_cast<const InternalData &> (fe_internal);
473 
474  // call FESystem's method to fill everything without enrichment function
475  fe_system->fill_fe_values(cell,
476  cell_similarity,
477  quadrature,
478  mapping,
479  mapping_internal,
480  mapping_data,
481  *fe_data.fesystem_data,
482  output_data);
483 
484  if (is_enriched)
485  multiply_by_enrichment(quadrature,
486  fe_data,
487  mapping_data,
488  cell,
489  output_data);
490 }
491 
492 
493 template <int dim, int spacedim>
494 void
497  const unsigned int face_no,
498  const Quadrature< dim-1 > &quadrature,
499  const Mapping< dim, spacedim > &mapping,
500  const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal,
501  const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data,
502  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
504 ) const
505 {
506  Assert (dynamic_cast<const InternalData *> (&fe_internal) != nullptr,
507  ExcInternalError());
508  const InternalData &fe_data = static_cast<const InternalData &> (fe_internal);
509 
510  // call FESystem's method to fill everything without enrichment function
511  fe_system->fill_fe_face_values(cell,
512  face_no,
513  quadrature,
514  mapping,
515  mapping_internal,
516  mapping_data,
517  *fe_data.fesystem_data,
518  output_data);
519 
520  if (is_enriched)
521  multiply_by_enrichment(quadrature,
522  fe_data,
523  mapping_data,
524  cell,
525  output_data);
526 }
527 
528 
529 template <int dim, int spacedim>
530 void
533  const unsigned int face_no,
534  const unsigned int sub_no,
535  const Quadrature< dim-1 > &quadrature,
536  const Mapping< dim, spacedim > &mapping,
537  const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal,
538  const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data,
539  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
541 ) const
542 {
543  Assert (dynamic_cast<const InternalData *> (&fe_internal) != nullptr,
544  ExcInternalError());
545  const InternalData &fe_data = static_cast<const InternalData &> (fe_internal);
546 
547  // call FESystem's method to fill everything without enrichment function
548  fe_system->fill_fe_subface_values(cell,
549  face_no,
550  sub_no,
551  quadrature,
552  mapping,
553  mapping_internal,
554  mapping_data,
555  *fe_data.fesystem_data,
556  output_data);
557 
558  if (is_enriched)
559  multiply_by_enrichment(quadrature,
560  fe_data,
561  mapping_data,
562  cell,
563  output_data);
564 
565 }
566 
567 
568 template <int dim, int spacedim>
569 template <int dim_1>
570 void
572 (const Quadrature<dim_1> &quadrature,
573  const InternalData &fe_data,
577 {
578  // mapping_data will contain quadrature points on the real element.
579  // fe_internal is needed to get update flags
580  // finally, output_data should store all the results we need.
581 
582  // Either dim_1==dim
583  // (fill_fe_values) or dim_1==dim-1
584  // (fill_fe_(sub)face_values)
585  Assert(dim_1==dim || dim_1==dim-1, ExcInternalError());
586  const UpdateFlags flags = fe_data.update_each;
587 
588  const unsigned int n_q_points = quadrature.size();
589 
590  // First, populate output_data object (that shall hold everything requested such
591  // as shape value, gradients, hessians, etc) from each base element.
592  // That is almost identical to FESystem::compute_fill_one_base(),
593  // the difference being that we do it irrespectively of cell_similarity
594  // and use base_fe_data.update_flags
595 
596  // TODO: do we need it only for dim_1 == dim (i.e. fill_fe_values)?
597  if (dim_1 == dim)
598  for (unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
599  {
601  base_fe = base_element(base_no);
603  base_fe_data = fe_data.get_fe_data(base_no);
605  base_data = fe_data.get_fe_output_object(base_no);
606 
607  const UpdateFlags base_flags = base_fe_data.update_each;
608 
609  for (unsigned int system_index=0; system_index<this->dofs_per_cell;
610  ++system_index)
611  if (this->system_to_base_table[system_index].first.first == base_no)
612  {
613  const unsigned int
614  base_index = this->system_to_base_table[system_index].second;
615  Assert (base_index<base_fe.dofs_per_cell, ExcInternalError());
616 
617  // now copy. if the shape function is primitive, then there
618  // is only one value to be copied, but for non-primitive
619  // elements, there might be more values to be copied
620  //
621  // so, find out from which index to take this one value, and
622  // to which index to put
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)
628  in_index += base_fe.n_nonzero_components(i);
629 
630  // then loop over the number of components to be copied
631  Assert (this->n_nonzero_components(system_index) ==
632  base_fe.n_nonzero_components(base_index),
633  ExcInternalError());
634  for (unsigned int s=0; s<this->n_nonzero_components(system_index); ++s)
635  {
636  if (base_flags & update_values)
637  for (unsigned int q=0; q<n_q_points; ++q)
638  output_data.shape_values[out_index+s][q] =
639  base_data.shape_values(in_index+s,q);
640 
641  if (base_flags & update_gradients)
642  for (unsigned int q=0; q<n_q_points; ++q)
643  output_data.shape_gradients[out_index+s][q] =
644  base_data.shape_gradients[in_index+s][q];
645 
646  if (base_flags & update_hessians)
647  for (unsigned int q=0; q<n_q_points; ++q)
648  output_data.shape_hessians[out_index+s][q] =
649  base_data.shape_hessians[in_index+s][q];
650  }
651  }
652  }
653 
654  Assert (base_no_mult_local_enriched_dofs.size() == fe_data.enrichment.size(),
655  ExcDimensionMismatch(base_no_mult_local_enriched_dofs.size(),
656  fe_data.enrichment.size()));
657  // calculate hessians, gradients and values for each function
658  for (unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
659  {
660  Assert (base_no_mult_local_enriched_dofs[base_no].size() == fe_data.enrichment[base_no].size(),
661  ExcDimensionMismatch(base_no_mult_local_enriched_dofs[base_no].size(),
662  fe_data.enrichment[base_no].size()));
663  for (unsigned int m=0; m < base_no_mult_local_enriched_dofs[base_no].size(); m++)
664  {
665  //Avoid evaluating quadrature points if no dofs are assigned. This
666  //happens when FE_Nothing is used together with other FE (i.e. FE_Q) as enrichments.
667  if (base_no_mult_local_enriched_dofs[base_no][m].size()==0)
668  continue;
669 
670  Assert (enrichments[base_no-1][m](cell) != nullptr,
671  ExcMessage("The pointer to the enrichment function is NULL"));
672 
673  Assert (enrichments[base_no-1][m](cell)->n_components == 1,
674  ExcMessage("Only scalar-valued enrichment functions are allowed"));
675 
676  if (flags & update_hessians)
677  {
678  Assert (fe_data.enrichment[base_no][m].hessians.size() == n_q_points,
679  ExcDimensionMismatch(fe_data.enrichment[base_no][m].hessians.size(),
680  n_q_points));
681  for (unsigned int q=0; q<n_q_points; q++)
682  fe_data.enrichment[base_no][m].hessians[q] = enrichments[base_no-1][m](cell)->hessian (mapping_data.quadrature_points[q]);
683  }
684 
685  if (flags & update_gradients)
686  {
687  Assert (fe_data.enrichment[base_no][m].gradients.size() == n_q_points,
688  ExcDimensionMismatch(fe_data.enrichment[base_no][m].gradients.size(),
689  n_q_points));
690  for (unsigned int q=0; q<n_q_points; q++)
691  fe_data.enrichment[base_no][m].gradients[q] = enrichments[base_no-1][m](cell)->gradient(mapping_data.quadrature_points[q]);
692  }
693 
694  if (flags & update_values)
695  {
696  Assert (fe_data.enrichment[base_no][m].values.size() == n_q_points,
697  ExcDimensionMismatch(fe_data.enrichment[base_no][m].values.size(),
698  n_q_points));
699  for (unsigned int q=0; q<n_q_points; q++)
700  fe_data.enrichment[base_no][m].values[q] = enrichments[base_no-1][m](cell)->value (mapping_data.quadrature_points[q]);
701  }
702  }
703  }
704 
705  // Finally, update the standard data stored in output_data
706  // by expanding the product rule for enrichment function.
707  // note that the order if important, namely
708  // output_data.shape_XYZ contains values of standard FEM and we overwrite
709  // it with the updated one in the following order: hessians -> gradients -> values
710  if (flags & update_hessians)
711  {
712  for (unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
713  {
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++)
716  {
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)
719  {
720  const Tensor<2, spacedim> grad_grad = outer_product(output_data.shape_gradients[enriched_dof][q],
721  fe_data.enrichment[base_no][m].gradients[q]);
722  const Tensor<2,spacedim,double> sym_grad_grad = symmetrize (grad_grad) * 2.0; // symmetrize does [s+s^T]/2
723 
724  output_data.shape_hessians [enriched_dof][q]*= fe_data.enrichment[base_no][m].values[q];
725  output_data.shape_hessians [enriched_dof][q]+=
726  sym_grad_grad +
727  output_data.shape_values [enriched_dof][q] * fe_data.enrichment[base_no][m].hessians[q];
728  }
729  }
730  }
731  }
732 
733  if (flags & update_gradients)
734  for (unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
735  {
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++)
738  {
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)
741  {
742  output_data.shape_gradients[enriched_dof][q]*= fe_data.enrichment[base_no][m].values[q];
743  output_data.shape_gradients[enriched_dof][q]+=
744  output_data.shape_values [enriched_dof][q] * fe_data.enrichment[base_no][m].gradients[q];
745  }
746  }
747  }
748 
749  if (flags & update_values)
750  for (unsigned int base_no = 1; base_no < this->n_base_elements(); base_no++)
751  {
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++)
754  {
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)
757  {
758  output_data.shape_values[enriched_dof][q] *= fe_data.enrichment[base_no][m].values[q];
759  }
760  }
761  }
762 }
763 
764 
765 template <int dim, int spacedim>
769 {
770  return *fe_system;
771 }
772 
773 
774 template <int dim, int spacedim>
775 bool
778 {
779  return true;
780 }
781 
782 
783 template <int dim, int spacedim>
784 void
787  FullMatrix<double> &matrix) const
788 {
789  if (const FE_Enriched<dim,spacedim> *fe_enr_other
790  = dynamic_cast<const FE_Enriched<dim,spacedim>*>(&source))
791  {
792  fe_system->get_face_interpolation_matrix(fe_enr_other->get_fe_system(),
793  matrix);
794  }
795  else
796  {
797  AssertThrow(false, (typename FiniteElement<dim,spacedim>::
798  ExcInterpolationNotImplemented()));
799  }
800 }
801 
802 
803 template <int dim, int spacedim>
804 void
807  const unsigned int subface,
808  FullMatrix<double> &matrix) const
809 {
810  if (const FE_Enriched<dim,spacedim> *fe_enr_other
811  = dynamic_cast<const FE_Enriched<dim,spacedim>*>(&source))
812  {
813  fe_system->get_subface_interpolation_matrix(fe_enr_other->get_fe_system(),
814  subface,
815  matrix);
816  }
817  else
818  {
819  AssertThrow(false,(typename FiniteElement<dim,spacedim>::
820  ExcInterpolationNotImplemented()));
821  }
822 }
823 
824 
825 template <int dim, int spacedim>
826 std::vector<std::pair<unsigned int, unsigned int> >
829 {
830  if (const FE_Enriched<dim,spacedim> *fe_enr_other
831  = dynamic_cast<const FE_Enriched<dim,spacedim>*>(&fe_other))
832  {
833  return fe_system->hp_vertex_dof_identities(fe_enr_other->get_fe_system());
834  }
835  else
836  {
837  Assert (false, ExcNotImplemented());
838  return std::vector<std::pair<unsigned int, unsigned int> >();
839  }
840 }
841 
842 
843 template <int dim, int spacedim>
844 std::vector<std::pair<unsigned int, unsigned int> >
847 {
848  if (const FE_Enriched<dim,spacedim> *fe_enr_other
849  = dynamic_cast<const FE_Enriched<dim,spacedim>*>(&fe_other))
850  {
851  return fe_system->hp_line_dof_identities(fe_enr_other->get_fe_system());
852  }
853  else
854  {
855  Assert (false, ExcNotImplemented());
856  return std::vector<std::pair<unsigned int, unsigned int> >();
857  }
858 }
859 
860 
861 template <int dim, int spacedim>
862 std::vector<std::pair<unsigned int, unsigned int> >
865 {
866  if (const FE_Enriched<dim,spacedim> *fe_enr_other
867  = dynamic_cast<const FE_Enriched<dim,spacedim>*>(&fe_other))
868  {
869  return fe_system->hp_quad_dof_identities(fe_enr_other->get_fe_system());
870  }
871  else
872  {
873  Assert (false, ExcNotImplemented());
874  return std::vector<std::pair<unsigned int, unsigned int> >();
875  }
876 }
877 
878 
879 template <int dim, int spacedim>
883 {
884  // need to decide which element constrain another.
885  // for example Q(2) dominate Q(4) and thus some DoFs of Q(4) will be constrained.
886  // If we have Q(2) and Q(4)+POU, then it's clear that Q(2) dominates,
887  // namely our DoFs will be constrained to make field continuous.
888  // However, we need to check for situations like Q(4) vs Q(2)+POU.
889  // In that case the domination for the underlying FEs should be the other way,
890  // but this implies that we can't constrain POU dofs to make the field continuous.
891  // In that case, through an error
892 
893  // if it's also enriched, do domination based on each one's FESystem
894  if (const FE_Enriched<dim,spacedim> *fe_enr_other
895  = dynamic_cast<const FE_Enriched<dim,spacedim>*>(&fe_other))
896  {
897  return fe_system->compare_for_face_domination(fe_enr_other->get_fe_system());
898  }
899  else
900  {
901  Assert (false, ExcNotImplemented());
903  }
904 }
905 
906 template <int dim, int spacedim>
907 const FullMatrix<double> &
909  const RefinementCase<dim> &refinement_case) const
910 {
911  return fe_system->get_prolongation_matrix(child, refinement_case);
912 }
913 
914 
915 template <int dim, int spacedim>
916 const FullMatrix<double> &
918  const RefinementCase<dim> &refinement_case) const
919 {
920  return fe_system->get_restriction_matrix(child, refinement_case);
921 }
922 
923 
924 /* ----------------------- FESystem::InternalData ------------------- */
925 
926 
927 template <int dim, int spacedim>
929  :
930  fesystem_data(std::move(fesystem_data))
931 {}
932 
933 
934 template <int dim, int spacedim>
937 InternalData::get_fe_data (const unsigned int base_no) const
938 {
939  return fesystem_data->get_fe_data(base_no);
940 }
941 
942 
943 template <int dim, int spacedim>
946 InternalData::get_fe_output_object (const unsigned int base_no) const
947 {
948  return fesystem_data->get_fe_output_object(base_no);
949 }
950 
951 
952 // explicit instantiations
953 #include "fe_enriched.inst"
954 
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
Definition: fe_enriched.cc:846
Shape function values.
bool is_dominating() const
Definition: fe_nothing.cc:176
Definition: tria.h:67
InternalData(std::unique_ptr< typename FESystem< dim, spacedim >::InternalData > fesystem_data)
Definition: fe_enriched.cc:928
unsigned int n_nonzero_components(const unsigned int i) const
Definition: fe.h:3217
FE_Enriched(const FiniteElement< dim, spacedim > &fe_base, const FiniteElement< dim, spacedim > &fe_enriched, const Function< spacedim > *enrichment_function)
Definition: fe_enriched.cc:133
void initialize(const std::vector< const FiniteElement< dim, spacedim > *> &fes, const std::vector< unsigned int > &multiplicities)
Definition: fe_enriched.cc:372
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
Definition: fe_enriched.cc:946
void build_face_tables(std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > &face_system_to_base_table, std::vector< std::pair< unsigned int, unsigned int > > &face_system_to_component_table, const FiniteElement< dim, spacedim > &finite_element, const bool do_tensor_product=true)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_enriched.cc:806
std::vector< std::vector< EnrichmentValues > > enrichment
Definition: fe_enriched.h:514
STL namespace.
Transformed quadrature points.
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe_enriched.cc:452
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
Definition: fe_enriched.cc:572
std::unique_ptr< To > dynamic_unique_cast(std::unique_ptr< From > &&p)
Definition: utilities.h:592
const FESystem< dim, spacedim > & get_fe_system() const
Definition: fe_enriched.cc:768
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
Definition: fe_enriched.cc:267
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual bool hp_constraints_are_implemented() const
Definition: fe_enriched.cc:777
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe_enriched.cc:238
Third derivatives of shape functions.
#define Assert(cond, exc)
Definition: exceptions.h:1142
UpdateFlags
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_enriched.cc:908
FiniteElement< dim, spacedim >::InternalDataBase & get_fe_data(const unsigned int base_no) const
Definition: fe_enriched.cc:937
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
virtual std::string get_name() const
Definition: fe_enriched.cc:429
Definition: fe.h:36
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)
Definition: utilities.cc:176
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
std::vector< Point< spacedim > > quadrature_points
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
Definition: fe_enriched.cc:331
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_enriched.cc:864
Shape function gradients.
std::vector< std::vector< std::vector< unsigned int > > > base_no_mult_local_enriched_dofs
Definition: fe_enriched.h:523
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_enriched.cc:828
void build_cell_tables(std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > &system_to_base_table, std::vector< std::pair< unsigned int, unsigned int > > &system_to_component_table, std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > &component_to_base_table, const FiniteElement< dim, spacedim > &finite_element, const bool do_tensor_product=true)
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_enriched.cc:786
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
Definition: fe_enriched.cc:345
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
Definition: fe_enriched.cc:292
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_enriched.cc:248
unsigned int n_base_elements() const
Definition: fe.h:3046
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_enriched.cc:917
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_enriched.cc:882
UpdateFlags update_each
Definition: fe.h:702
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
Definition: fe.h:2499
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
Definition: fe_enriched.cc:359