59template <
int spacedim>
60template <
typename InputVector>
69 const InputVector & solution,
73 const unsigned int n_threads,
79 const std::vector<const InputVector *> solutions(1, &solution);
80 std::vector<Vector<float> *> errors(1, &error);
97template <
int spacedim>
98template <
typename InputVector>
106 const InputVector & solution,
110 const unsigned int n_threads,
131template <
int spacedim>
132template <
typename InputVector>
140 const std::vector<const InputVector *> &solutions,
144 const unsigned int n_threads,
165template <
int spacedim>
166template <
typename InputVector>
175 const InputVector & solution,
179 const unsigned int n_threads,
185 const std::vector<const InputVector *> solutions(1, &solution);
186 std::vector<Vector<float> *> errors(1, &error);
202template <
int spacedim>
203template <
typename InputVector>
211 const InputVector & solution,
215 const unsigned int n_threads,
236template <
int spacedim>
237template <
typename InputVector>
245 const std::vector<const InputVector *> &solutions,
249 const unsigned int n_threads,
270template <
int spacedim>
271template <
typename InputVector>
280 const std::vector<const InputVector *> & ,
294template <
int spacedim>
295template <
typename InputVector>
304 const std::vector<const InputVector *> &solutions,
314 using number =
typename InputVector::value_type;
315#ifdef DEAL_II_WITH_P4EST
325 "For parallel distributed triangulations, the only "
326 "valid subdomain_id that can be passed here is the "
327 "one that corresponds to the locally owned subdomain id."));
334 .locally_owned_subdomain() :
341 const unsigned int n_solution_vectors = solutions.size();
346 ExcMessage(
"You are not allowed to list the special boundary "
347 "indicator for internal boundaries in your boundary "
350 for (
const auto &boundary_function : neumann_bc)
352 (void)boundary_function;
353 Assert(boundary_function.second->n_components == n_components,
355 boundary_function.second->n_components,
364 Assert((coefficient ==
nullptr) ||
370 Assert(solutions.size() == errors.size(),
372 for (
unsigned int n = 0; n < solutions.size(); ++n)
376 Assert((coefficient ==
nullptr) ||
381 for (
const auto &boundary_function : neumann_bc)
383 (void)boundary_function;
384 Assert(boundary_function.second->n_components == n_components,
386 boundary_function.second->n_components,
391 for (
unsigned int n = 0; n < n_solution_vectors; ++n)
398 std::vector<std::vector<std::vector<Tensor<1, spacedim, number>>>>
399 gradients_here(n_solution_vectors,
403 std::vector<std::vector<std::vector<Tensor<1, spacedim, number>>>>
404 gradients_neighbor(gradients_here);
405 std::vector<Vector<typename ProductType<number, double>::type>>
406 grad_dot_n_neighbor(n_solution_vectors,
414 if (coefficient ==
nullptr)
415 for (
unsigned int c = 0; c < n_components; ++c)
416 coefficient_values(c) = 1;
444 for (
unsigned int n = 0; n < n_solution_vectors; ++n)
445 (*errors[n])(cell->active_cell_index()) = 0;
448 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
450 *solutions[s], gradients_here[s]);
454 for (
unsigned int n = 0; n < 2; ++n)
457 auto neighbor = cell->neighbor(n);
459 while (neighbor->has_children())
460 neighbor = neighbor->child(n == 0 ? 1 : 0);
462 fe_face_values.
reinit(cell, n);
468 fe_values.
reinit(neighbor);
470 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
472 *solutions[s], gradients_neighbor[s]);
474 fe_face_values.
reinit(neighbor, n == 0 ? 1 : 0);
481 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
482 for (
unsigned int c = 0; c < n_components; ++c)
483 grad_dot_n_neighbor[s](c) =
484 -(gradients_neighbor[s][n == 0 ? 1 : 0][c] *
487 else if (neumann_bc.find(n) != neumann_bc.end())
491 if (n_components == 1)
493 const typename InputVector::value_type v =
494 neumann_bc.find(n)->second->value(cell->vertex(n));
496 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
497 grad_dot_n_neighbor[s](0) = v;
502 neumann_bc.find(n)->second->vector_value(cell->vertex(n),
505 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
506 grad_dot_n_neighbor[s] = v;
511 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
512 grad_dot_n_neighbor[s] = 0;
516 if (coefficient !=
nullptr)
520 const double c_value = coefficient->
value(cell->vertex(n));
521 for (
unsigned int c = 0; c < n_components; ++c)
522 coefficient_values(c) = c_value;
530 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
531 for (
unsigned int component = 0; component < n_components;
533 if (component_mask[component] ==
true)
538 gradients_here[s][n][component] * normal;
541 ((grad_dot_n_here - grad_dot_n_neighbor[s](component)) *
542 coefficient_values(component));
543 (*errors[s])(cell->active_cell_index()) +=
546 double>::type>::abs_square(jump) *
551 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
552 (*errors[s])(cell->active_cell_index()) =
553 std::sqrt((*errors[s])(cell->active_cell_index()));
559#include "error_estimator_1d.inst"
bool represents_n_components(const unsigned int n) const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
unsigned int n_components() const
const unsigned int n_components
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
@ cell_diameter_over_24
Kelly error estimator with the factor .
Abstract base class for mapping classes.
unsigned int n_active_cells() const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell, const unsigned int face_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
const FEValuesType & get_present_fe_values() const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
void push_back(const Mapping< dim, spacedim > &new_mapping)
types::subdomain_id locally_owned_subdomain() const override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_normal_vectors
Normal vectors.
@ update_gradients
Shape function gradients.
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcInvalidComponentMask()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcInvalidCoefficient()
static ::ExceptionBase & ExcNoSolutions()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ valid
Iterator points to a valid object.
const types::material_id invalid_material_id
const types::boundary_id internal_face_boundary_id
const types::subdomain_id invalid_subdomain_id
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
unsigned int subdomain_id
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type