16 #include <deal.II/base/geometry_info.h>
17 #include <deal.II/base/quadrature.h>
18 #include <deal.II/base/quadrature_lib.h>
19 #include <deal.II/base/thread_management.h>
20 #include <deal.II/base/work_stream.h>
22 #include <deal.II/distributed/tria.h>
24 #include <deal.II/dofs/dof_accessor.h>
25 #include <deal.II/dofs/dof_handler.h>
27 #include <deal.II/fe/fe.h>
28 #include <deal.II/fe/fe_update_flags.h>
29 #include <deal.II/fe/fe_values.h>
30 #include <deal.II/fe/mapping_q1.h>
32 #include <deal.II/grid/tria_iterator.h>
34 #include <deal.II/hp/fe_values.h>
35 #include <deal.II/hp/mapping_collection.h>
36 #include <deal.II/hp/q_collection.h>
38 #include <deal.II/lac/block_vector.h>
39 #include <deal.II/lac/la_parallel_block_vector.h>
40 #include <deal.II/lac/la_parallel_vector.h>
41 #include <deal.II/lac/la_vector.h>
42 #include <deal.II/lac/petsc_block_vector.h>
43 #include <deal.II/lac/petsc_vector.h>
44 #include <deal.II/lac/trilinos_parallel_block_vector.h>
45 #include <deal.II/lac/trilinos_vector.h>
46 #include <deal.II/lac/vector.h>
48 #include <deal.II/numerics/error_estimator.h>
56 DEAL_II_NAMESPACE_OPEN
60 template <
int spacedim>
61 template <
typename InputVector,
typename DoFHandlerType>
65 const DoFHandlerType & dof_handler,
70 const InputVector & solution,
71 Vector<float> & error,
74 const unsigned int n_threads,
80 const std::vector<const InputVector *> solutions(1, &solution);
81 std::vector<Vector<float> *> errors(1, &error);
98 template <
int spacedim>
99 template <
typename InputVector,
typename DoFHandlerType>
102 const DoFHandlerType &dof_handler,
107 const InputVector & solution,
108 Vector<float> & error,
111 const unsigned int n_threads,
132 template <
int spacedim>
133 template <
typename InputVector,
typename DoFHandlerType>
136 const DoFHandlerType &dof_handler,
141 const std::vector<const InputVector *> &solutions,
142 std::vector<Vector<float> *> & errors,
145 const unsigned int n_threads,
166 template <
int spacedim>
167 template <
typename InputVector,
typename DoFHandlerType>
171 const DoFHandlerType & dof_handler,
176 const InputVector & solution,
177 Vector<float> & error,
180 const unsigned int n_threads,
186 const std::vector<const InputVector *> solutions(1, &solution);
187 std::vector<Vector<float> *> errors(1, &error);
203 template <
int spacedim>
204 template <
typename InputVector,
typename DoFHandlerType>
207 const DoFHandlerType & dof_handler,
212 const InputVector & solution,
213 Vector<float> & error,
216 const unsigned int n_threads,
237 template <
int spacedim>
238 template <
typename InputVector,
typename DoFHandlerType>
241 const DoFHandlerType & dof_handler,
246 const std::vector<const InputVector *> &solutions,
247 std::vector<Vector<float> *> & errors,
250 const unsigned int n_threads,
271 template <
int spacedim>
272 template <
typename InputVector,
typename DoFHandlerType>
276 const DoFHandlerType & ,
281 const std::vector<const InputVector *> & ,
282 std::vector<Vector<float> *> & ,
295 template <
int spacedim>
296 template <
typename InputVector,
typename DoFHandlerType>
300 const DoFHandlerType & dof_handler,
305 const std::vector<const InputVector *> &solutions,
306 std::vector<Vector<float> *> & errors,
315 using number =
typename InputVector::value_type;
316 #ifdef DEAL_II_WITH_P4EST
318 &dof_handler.get_triangulation()) !=
nullptr)
323 dof_handler.get_triangulation())
324 .locally_owned_subdomain()),
326 "For parallel distributed triangulations, the only "
327 "valid subdomain_id that can be passed here is the "
328 "one that corresponds to the locally owned subdomain id."));
332 &dof_handler.get_triangulation()) !=
nullptr) ?
334 dof_handler.get_triangulation())
335 .locally_owned_subdomain() :
341 const unsigned int n_components = dof_handler.get_fe(0).n_components();
342 const unsigned int n_solution_vectors = solutions.size();
347 ExcMessage(
"You are not allowed to list the special boundary "
348 "indicator for internal boundaries in your boundary "
351 for (
const auto &boundary_function : neumann_bc)
353 (void)boundary_function;
354 Assert(boundary_function.second->n_components == n_components,
356 boundary_function.second->n_components,
365 Assert((coefficient ==
nullptr) ||
371 Assert(solutions.size() == errors.size(),
373 for (
unsigned int n = 0; n < solutions.size(); ++n)
374 Assert(solutions[n]->size() == dof_handler.n_dofs(),
377 Assert((coefficient ==
nullptr) ||
382 for (
const auto &boundary_function : neumann_bc)
384 (void)boundary_function;
385 Assert(boundary_function.second->n_components == n_components,
387 boundary_function.second->n_components,
392 for (
unsigned int n = 0; n < n_solution_vectors; ++n)
393 (*errors[n]).reinit(dof_handler.get_triangulation().n_active_cells());
399 std::vector<std::vector<std::vector<Tensor<1, spacedim, number>>>>
400 gradients_here(n_solution_vectors,
404 std::vector<std::vector<std::vector<Tensor<1, spacedim, number>>>>
405 gradients_neighbor(gradients_here);
406 std::vector<Vector<typename ProductType<number, double>::type>>
407 grad_dot_n_neighbor(n_solution_vectors,
408 Vector<
typename ProductType<number, double>::type>(
415 if (coefficient ==
nullptr)
416 for (
unsigned int c = 0; c < n_components; ++c)
417 coefficient_values(c) = 1;
434 fe, q_face_collection, update_normal_vectors);
439 for (
typename DoFHandlerType::active_cell_iterator cell =
440 dof_handler.begin_active();
441 cell != dof_handler.end();
444 (cell->subdomain_id() == subdomain_id)) &&
446 (cell->material_id() == material_id)))
448 for (
unsigned int n = 0; n < n_solution_vectors; ++n)
449 (*errors[n])(cell->active_cell_index()) = 0;
452 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
453 fe_values.get_present_fe_values().get_function_gradients(
454 *solutions[s], gradients_here[s]);
458 for (
unsigned int n = 0; n < 2; ++n)
461 typename DoFHandlerType::cell_iterator neighbor = cell->neighbor(n);
463 while (neighbor->has_children())
464 neighbor = neighbor->child(n == 0 ? 1 : 0);
466 fe_face_values.
reinit(cell, n);
468 .get_all_normal_vectors()[0];
472 fe_values.
reinit(neighbor);
474 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
475 fe_values.get_present_fe_values().get_function_gradients(
476 *solutions[s], gradients_neighbor[s]);
478 fe_face_values.
reinit(neighbor, n == 0 ? 1 : 0);
480 fe_face_values.get_present_fe_values()
481 .get_all_normal_vectors()[0];
485 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
486 for (
unsigned int c = 0; c < n_components; ++c)
487 grad_dot_n_neighbor[s](c) =
488 -(gradients_neighbor[s][n == 0 ? 1 : 0][c] *
491 else if (neumann_bc.find(n) != neumann_bc.end())
495 if (n_components == 1)
497 const typename InputVector::value_type v =
498 neumann_bc.find(n)->second->value(cell->vertex(n));
500 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
501 grad_dot_n_neighbor[s](0) = v;
505 Vector<typename InputVector::value_type> v(n_components);
506 neumann_bc.find(n)->second->vector_value(cell->vertex(n),
509 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
510 grad_dot_n_neighbor[s] = v;
515 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
516 grad_dot_n_neighbor[s] = 0;
520 if (coefficient !=
nullptr)
524 const double c_value = coefficient->
value(cell->vertex(n));
525 for (
unsigned int c = 0; c < n_components; ++c)
526 coefficient_values(c) = c_value;
534 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
535 for (
unsigned int component = 0; component < n_components;
537 if (component_mask[component] ==
true)
540 const typename ProductType<number, double>::type
542 gradients_here[s][n][component] * normal;
544 const typename ProductType<number, double>::type jump =
545 ((grad_dot_n_here - grad_dot_n_neighbor[s](component)) *
546 coefficient_values(component));
547 (*errors[s])(cell->active_cell_index()) +=
550 double>::type>::abs_square(jump) *
555 for (
unsigned int s = 0; s < n_solution_vectors; ++s)
556 (*errors[s])(cell->active_cell_index()) =
557 std::sqrt((*errors[s])(cell->active_cell_index()));
563 #include "error_estimator_1d.inst"
566 DEAL_II_NAMESPACE_CLOSE