deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
error_estimator_1d.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2013 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
19
21
24
25#include <deal.II/fe/fe.h>
28#include <deal.II/fe/mapping.h>
29
31
35
43#include <deal.II/lac/vector.h>
44
46
47#include <algorithm>
48#include <cmath>
49#include <functional>
50#include <numeric>
51#include <vector>
52
54
55
56template <int spacedim>
57template <typename Number>
58void
61 const DoFHandler<1, spacedim> &dof_handler,
62 const Quadrature<0> &quadrature,
63 const std::map<types::boundary_id, const Function<spacedim, Number> *>
64 &neumann_bc,
65 const ReadVector<Number> &solution,
66 Vector<float> &error,
67 const ComponentMask &component_mask,
68 const Function<spacedim> *coefficients,
69 const unsigned int n_threads,
70 const types::subdomain_id subdomain_id,
71 const types::material_id material_id,
72 const Strategy strategy)
73{
74 // just pass on to the other function
75 std::vector<const ReadVector<Number> *> solutions(1, &solution);
76 std::vector<Vector<float> *> errors(1, &error);
77 ArrayView<Vector<float> *> error_view = make_array_view(errors);
79 dof_handler,
80 quadrature,
81 neumann_bc,
82 make_array_view(solutions),
83 error_view,
84 component_mask,
85 coefficients,
86 n_threads,
87 subdomain_id,
88 material_id,
89 strategy);
90}
91
92
93
94template <int spacedim>
95template <typename Number>
96void
98 const DoFHandler<1, spacedim> &dof_handler,
99 const Quadrature<0> &quadrature,
100 const std::map<types::boundary_id, const Function<spacedim, Number> *>
101 &neumann_bc,
102 const ReadVector<Number> &solution,
103 Vector<float> &error,
104 const ComponentMask &component_mask,
105 const Function<spacedim> *coefficients,
106 const unsigned int n_threads,
107 const types::subdomain_id subdomain_id,
108 const types::material_id material_id,
109 const Strategy strategy)
110{
111 const auto reference_cell = ReferenceCells::Line;
112 estimate(reference_cell.template get_default_linear_mapping<1, spacedim>(),
113 dof_handler,
114 quadrature,
115 neumann_bc,
116 solution,
117 error,
118 component_mask,
119 coefficients,
120 n_threads,
121 subdomain_id,
122 material_id,
123 strategy);
124}
125
126
127
128template <int spacedim>
129template <typename Number>
130void
132 const DoFHandler<1, spacedim> &dof_handler,
133 const Quadrature<0> &quadrature,
134 const std::map<types::boundary_id, const Function<spacedim, Number> *>
135 &neumann_bc,
136 const ArrayView<const ReadVector<Number> *> &solutions,
137 ArrayView<Vector<float> *> &errors,
138 const ComponentMask &component_mask,
139 const Function<spacedim> *coefficients,
140 const unsigned int n_threads,
141 const types::subdomain_id subdomain_id,
142 const types::material_id material_id,
143 const Strategy strategy)
144{
145 const auto reference_cell = ReferenceCells::Line;
146 estimate(reference_cell.template get_default_linear_mapping<1, spacedim>(),
147 dof_handler,
148 quadrature,
149 neumann_bc,
150 solutions,
151 errors,
152 component_mask,
153 coefficients,
154 n_threads,
155 subdomain_id,
156 material_id,
157 strategy);
158}
159
160
161
162template <int spacedim>
163template <typename Number>
164void
167 const DoFHandler<1, spacedim> &dof_handler,
168 const hp::QCollection<0> &quadrature,
169 const std::map<types::boundary_id, const Function<spacedim, Number> *>
170 &neumann_bc,
171 const ReadVector<Number> &solution,
172 Vector<float> &error,
173 const ComponentMask &component_mask,
174 const Function<spacedim> *coefficients,
175 const unsigned int n_threads,
176 const types::subdomain_id subdomain_id,
177 const types::material_id material_id,
178 const Strategy strategy)
179{
180 // just pass on to the other function
181 std::vector<const ReadVector<Number> *> solutions(1, &solution);
182 std::vector<Vector<float> *> errors(1, &error);
183 ArrayView<Vector<float> *> error_view = make_array_view(errors);
185 dof_handler,
186 quadrature,
187 neumann_bc,
188 make_array_view(solutions),
189 error_view,
190 component_mask,
191 coefficients,
192 n_threads,
193 subdomain_id,
194 material_id,
195 strategy);
196}
197
198
199
200template <int spacedim>
201template <typename Number>
202void
205 const DoFHandler<1, spacedim> &dof_handler,
206 const hp::QCollection<0> &quadrature,
207 const std::map<types::boundary_id, const Function<spacedim, Number> *>
208 &neumann_bc,
209 const ReadVector<Number> &solution,
210 Vector<float> &error,
211 const ComponentMask &component_mask,
212 const Function<spacedim> *coefficients,
213 const unsigned int n_threads,
214 const types::subdomain_id subdomain_id,
215 const types::material_id material_id,
216 const Strategy strategy)
217{
218 // DEPRECATED
219 // just pass on to the other function
220 std::vector<const ReadVector<Number> *> solutions(1, &solution);
221 std::vector<Vector<float> *> errors(1, &error);
222 ArrayView<Vector<float> *> error_view = make_array_view(errors);
223 const hp::MappingCollection<1, spacedim> mapping_collection(mapping);
224 estimate(mapping_collection,
225 dof_handler,
226 quadrature,
227 neumann_bc,
228 make_array_view(solutions),
229 error_view,
230 component_mask,
231 coefficients,
232 n_threads,
233 subdomain_id,
234 material_id,
235 strategy);
236}
237
238
239
240template <int spacedim>
241template <typename Number>
242void
244 const DoFHandler<1, spacedim> &dof_handler,
245 const hp::QCollection<0> &quadrature,
246 const std::map<types::boundary_id, const Function<spacedim, Number> *>
247 &neumann_bc,
248 const ReadVector<Number> &solution,
249 Vector<float> &error,
250 const ComponentMask &component_mask,
251 const Function<spacedim> *coefficients,
252 const unsigned int n_threads,
253 const types::subdomain_id subdomain_id,
254 const types::material_id material_id,
255 const Strategy strategy)
256{
257 const auto reference_cell = ReferenceCells::Line;
259 reference_cell.template get_default_linear_mapping<1, spacedim>());
261 dof_handler,
262 quadrature,
263 neumann_bc,
264 solution,
265 error,
266 component_mask,
267 coefficients,
268 n_threads,
269 subdomain_id,
270 material_id,
271 strategy);
272}
273
274
275
276template <int spacedim>
277template <typename Number>
278void
280 const DoFHandler<1, spacedim> &dof_handler,
281 const hp::QCollection<0> &quadrature,
282 const std::map<types::boundary_id, const Function<spacedim, Number> *>
283 &neumann_bc,
284 const ArrayView<const ReadVector<Number> *> &solutions,
285 ArrayView<Vector<float> *> &errors,
286 const ComponentMask &component_mask,
287 const Function<spacedim> *coefficients,
288 const unsigned int n_threads,
289 const types::subdomain_id subdomain_id,
290 const types::material_id material_id,
291 const Strategy strategy)
292{
293 const auto reference_cell = ReferenceCells::Line;
295 reference_cell.template get_default_linear_mapping<1, spacedim>());
297 dof_handler,
298 quadrature,
299 neumann_bc,
300 solutions,
301 errors,
302 component_mask,
303 coefficients,
304 n_threads,
305 subdomain_id,
306 material_id,
307 strategy);
308}
309
310
311
312template <int spacedim>
313template <typename Number>
314void
317 const DoFHandler<1, spacedim> &dof_handler,
318 const hp::QCollection<0> &,
319 const std::map<types::boundary_id, const Function<spacedim, Number> *>
320 &neumann_bc,
321 const ArrayView<const ReadVector<Number> *> &solutions,
322 ArrayView<Vector<float> *> &errors,
323 const ComponentMask &component_mask,
324 const Function<spacedim> *coefficient,
325 const unsigned int,
326 const types::subdomain_id subdomain_id_,
327 const types::material_id material_id,
328 const Strategy strategy)
329{
331 using number = Number;
333 if (const auto *triangulation = dynamic_cast<
335 &dof_handler.get_triangulation()))
336 {
337 Assert((subdomain_id_ == numbers::invalid_subdomain_id) ||
338 (subdomain_id_ == triangulation->locally_owned_subdomain()),
340 "For distributed Triangulation objects and associated "
341 "DoFHandler objects, asking for any subdomain other than the "
342 "locally owned one does not make sense."));
343 subdomain_id = triangulation->locally_owned_subdomain();
344 }
345 else
346 {
347 subdomain_id = subdomain_id_;
348 }
349
350 const unsigned int n_components = dof_handler.get_fe(0).n_components();
351 const unsigned int n_solution_vectors = solutions.size();
352
353 // sanity checks
355 neumann_bc.end(),
356 ExcMessage("You are not allowed to list the special boundary "
357 "indicator for internal boundaries in your boundary "
358 "value map."));
359
360 for (const auto &boundary_function : neumann_bc)
361 {
362 (void)boundary_function;
363 Assert(boundary_function.second->n_components == n_components,
364 ExcInvalidBoundaryFunction(boundary_function.first,
365 boundary_function.second->n_components,
366 n_components));
367 }
368
369 Assert(component_mask.represents_n_components(n_components),
371 Assert(component_mask.n_selected_components(n_components) > 0,
373
374 Assert((coefficient == nullptr) ||
375 (coefficient->n_components == n_components) ||
376 (coefficient->n_components == 1),
378
379 Assert(solutions.size() > 0, ExcNoSolutions());
380 Assert(solutions.size() == errors.size(),
381 ExcIncompatibleNumberOfElements(solutions.size(), errors.size()));
382 for (unsigned int n = 0; n < solutions.size(); ++n)
383 Assert(solutions[n]->size() == dof_handler.n_dofs(),
384 ExcDimensionMismatch(solutions[n]->size(), dof_handler.n_dofs()));
385
386 Assert((coefficient == nullptr) ||
387 (coefficient->n_components == n_components) ||
388 (coefficient->n_components == 1),
390
391 for (const auto &boundary_function : neumann_bc)
392 {
393 (void)boundary_function;
394 Assert(boundary_function.second->n_components == n_components,
395 ExcInvalidBoundaryFunction(boundary_function.first,
396 boundary_function.second->n_components,
397 n_components));
398 }
399
400 // reserve one slot for each cell and set it to zero
401 for (unsigned int n = 0; n < n_solution_vectors; ++n)
402 (*errors[n]).reinit(dof_handler.get_triangulation().n_active_cells());
403
404 // fields to get the gradients on the present and the neighbor cell.
405 //
406 // for the neighbor gradient, we need several auxiliary fields, depending on
407 // the way we get it (see below)
408 std::vector<std::vector<std::vector<Tensor<1, spacedim, number>>>>
409 gradients_here(n_solution_vectors,
410 std::vector<std::vector<Tensor<1, spacedim, number>>>(
411 2,
412 std::vector<Tensor<1, spacedim, number>>(n_components)));
413 std::vector<std::vector<std::vector<Tensor<1, spacedim, number>>>>
414 gradients_neighbor(gradients_here);
415 std::vector<Vector<typename ProductType<number, double>::type>>
416 grad_dot_n_neighbor(n_solution_vectors,
418 n_components));
419
420 // reserve some space for coefficient values at one point. if there is no
421 // coefficient, then we fill it by unity once and for all and don't set it
422 // any more
423 Vector<double> coefficient_values(n_components);
424 if (coefficient == nullptr)
425 for (unsigned int c = 0; c < n_components; ++c)
426 coefficient_values(c) = 1;
427
428 const QTrapezoid<1> quadrature;
429 const hp::QCollection<1> q_collection(quadrature);
430 const QGauss<0> face_quadrature(1);
431 const hp::QCollection<0> q_face_collection(face_quadrature);
432
433 const hp::FECollection<1, spacedim> &fe = dof_handler.get_fe_collection();
434
435
437 fe,
438 q_collection,
440 hp::FEFaceValues<1, spacedim> fe_face_values(
441 /*mapping,*/ fe, q_face_collection, update_normal_vectors);
442
443 // loop over all cells and do something on the cells which we're told to
444 // work on. note that the error indicator is only a sum over the two
445 // contributions from the two vertices of each cell.
446 for (const auto &cell : dof_handler.active_cell_iterators())
447 if (((subdomain_id == numbers::invalid_subdomain_id) ||
448 (cell->subdomain_id() == subdomain_id)) &&
449 ((material_id == numbers::invalid_material_id) ||
450 (cell->material_id() == material_id)))
451 {
452 for (unsigned int n = 0; n < n_solution_vectors; ++n)
453 (*errors[n])(cell->active_cell_index()) = 0;
454
455 fe_values.reinit(cell);
456 for (unsigned int s = 0; s < n_solution_vectors; ++s)
457 fe_values.get_present_fe_values().get_function_gradients(
458 *solutions[s], gradients_here[s]);
459
460 // loop over the two points bounding this line. n==0 is left point,
461 // n==1 is right point
462 for (unsigned int n = 0; n < 2; ++n)
463 {
464 // find left or right active neighbor
465 auto neighbor = cell->neighbor(n);
466 if (neighbor.state() == IteratorState::valid)
467 while (neighbor->has_children())
468 neighbor = neighbor->child(n == 0 ? 1 : 0);
469
470 fe_face_values.reinit(cell, n);
471 Tensor<1, spacedim> normal =
472 fe_face_values.get_present_fe_values().get_normal_vectors()[0];
473
474 if (neighbor.state() == IteratorState::valid)
475 {
476 fe_values.reinit(neighbor);
477
478 for (unsigned int s = 0; s < n_solution_vectors; ++s)
479 fe_values.get_present_fe_values().get_function_gradients(
480 *solutions[s], gradients_neighbor[s]);
481
482 fe_face_values.reinit(neighbor, n == 0 ? 1 : 0);
483 Tensor<1, spacedim> neighbor_normal =
484 fe_face_values.get_present_fe_values()
485 .get_normal_vectors()[0];
486
487 // extract the gradient in normal direction of all the
488 // components.
489 for (unsigned int s = 0; s < n_solution_vectors; ++s)
490 for (unsigned int c = 0; c < n_components; ++c)
491 grad_dot_n_neighbor[s](c) =
492 -(gradients_neighbor[s][n == 0 ? 1 : 0][c] *
493 neighbor_normal);
494 }
495 else if (neumann_bc.find(n) != neumann_bc.end())
496 // if Neumann b.c., then fill the gradients field which will be
497 // used later on.
498 {
499 if (n_components == 1)
500 {
501 const Number v =
502 neumann_bc.find(n)->second->value(cell->vertex(n));
503
504 for (unsigned int s = 0; s < n_solution_vectors; ++s)
505 grad_dot_n_neighbor[s](0) = v;
506 }
507 else
508 {
509 Vector<Number> v(n_components);
510 neumann_bc.find(n)->second->vector_value(cell->vertex(n),
511 v);
512
513 for (unsigned int s = 0; s < n_solution_vectors; ++s)
514 grad_dot_n_neighbor[s] = v;
515 }
516 }
517 else
518 // fill with zeroes.
519 for (unsigned int s = 0; s < n_solution_vectors; ++s)
520 grad_dot_n_neighbor[s] = 0;
521
522 // if there is a coefficient, then evaluate it at the present
523 // position. if there is none, reuse the preset values.
524 if (coefficient != nullptr)
525 {
526 if (coefficient->n_components == 1)
527 {
528 const double c_value = coefficient->value(cell->vertex(n));
529 for (unsigned int c = 0; c < n_components; ++c)
530 coefficient_values(c) = c_value;
531 }
532 else
533 coefficient->vector_value(cell->vertex(n),
534 coefficient_values);
535 }
536
537
538 for (unsigned int s = 0; s < n_solution_vectors; ++s)
539 for (unsigned int component = 0; component < n_components;
540 ++component)
541 if (component_mask[component] == true)
542 {
543 // get gradient here
545 grad_dot_n_here =
546 gradients_here[s][n][component] * normal;
547
548 const typename ProductType<number, double>::type jump =
549 ((grad_dot_n_here - grad_dot_n_neighbor[s](component)) *
550 coefficient_values(component));
551 (*errors[s])(cell->active_cell_index()) +=
553 typename ProductType<number,
554 double>::type>::abs_square(jump) *
555 cell->diameter();
556 }
557 }
558
559 for (unsigned int s = 0; s < n_solution_vectors; ++s)
560 (*errors[s])(cell->active_cell_index()) =
561 std::sqrt((*errors[s])(cell->active_cell_index()));
562 }
563}
564
565
566
567template <int spacedim>
568template <typename Number>
569void
572 const DoFHandler<1, spacedim> &dof_handler,
573 const hp::QCollection<0> &quadrature,
574 const std::map<types::boundary_id, const Function<spacedim, Number> *>
575 &neumann_bc,
576 const ArrayView<const ReadVector<Number> *> &solutions,
577 ArrayView<Vector<float> *> &errors,
578 const ComponentMask &component_mask,
579 const Function<spacedim> *coefficients,
580 const unsigned int n_threads,
581 const types::subdomain_id subdomain_id,
582 const types::material_id material_id,
583 const Strategy strategy)
584{
585 // DEPRECATED
586 // just pass on to the other function
587 const hp::MappingCollection<1, spacedim> mapping_collection(mapping);
588 estimate(mapping_collection,
589 dof_handler,
590 quadrature,
591 neumann_bc,
592 solutions,
593 errors,
594 component_mask,
595 coefficients,
596 n_threads,
597 subdomain_id,
598 material_id,
599 strategy);
600}
601
602
603
604template <int spacedim>
605template <typename Number>
606void
609 const DoFHandler<1, spacedim> &dof_handler,
610 const Quadrature<0> &quadrature,
611 const std::map<types::boundary_id, const Function<spacedim, Number> *>
612 &neumann_bc,
613 const ArrayView<const ReadVector<Number> *> &solutions,
614 ArrayView<Vector<float> *> &errors,
615 const ComponentMask &component_mask,
616 const Function<spacedim> *coefficients,
617 const unsigned int n_threads,
618 const types::subdomain_id subdomain_id,
619 const types::material_id material_id,
620 const Strategy strategy)
621{
622 const hp::MappingCollection<1, spacedim> mapping_collection(mapping);
623 const hp::QCollection<0> quadrature_collection(quadrature);
624 estimate(mapping_collection,
625 dof_handler,
626 quadrature_collection,
627 neumann_bc,
628 solutions,
629 errors,
630 component_mask,
631 coefficients,
632 n_threads,
633 subdomain_id,
634 material_id,
635 strategy);
636}
637
638
639
640// explicit instantiations
641#include "numerics/error_estimator_1d.inst"
642
643
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:954
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 hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
unsigned int n_components() const
const unsigned int n_components
Definition function.h:164
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
@ cell_diameter_over_24
Kelly error estimator with the factor .
static void estimate(const Mapping< 1, spacedim > &mapping, const DoFHandler< 1, spacedim > &dof, const Quadrature< 0 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, const Function< spacedim > *coefficient=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)
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, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, 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)
Abstract base class for mapping classes.
Definition mapping.h:320
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)
Definition fe_values.cc:426
const FEValuesType & get_present_fe_values() const
Definition fe_values.h:695
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)
Definition fe_values.cc:296
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcNoSolutions()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcInvalidCoefficient()
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcInvalidComponentMask()
#define AssertThrow(cond, exc)
@ update_normal_vectors
Normal vectors.
@ update_gradients
Shape function gradients.
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:316
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
@ valid
Iterator points to a valid object.
constexpr ReferenceCell Line
constexpr types::boundary_id internal_face_boundary_id
Definition types.h:329
constexpr types::material_id invalid_material_id
Definition types.h:294
constexpr types::subdomain_id invalid_subdomain_id
Definition types.h:381
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
unsigned int material_id
Definition types.h:184
unsigned int boundary_id
Definition types.h:161
unsigned int subdomain_id
Definition types.h:52
typename internal::ProductTypeImpl< std::decay_t< T >, std::decay_t< U > >::type type