Reference documentation for deal.II version 8.4.2
error_estimator_1d.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2016 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 #include <deal.II/base/thread_management.h>
17 #include <deal.II/base/quadrature.h>
18 #include <deal.II/base/quadrature_lib.h>
19 #include <deal.II/base/work_stream.h>
20 #include <deal.II/lac/vector.h>
21 #include <deal.II/lac/parallel_vector.h>
22 #include <deal.II/lac/block_vector.h>
23 #include <deal.II/lac/parallel_block_vector.h>
24 #include <deal.II/lac/petsc_vector.h>
25 #include <deal.II/lac/petsc_block_vector.h>
26 #include <deal.II/lac/trilinos_vector.h>
27 #include <deal.II/lac/trilinos_block_vector.h>
28 #include <deal.II/grid/tria_iterator.h>
29 #include <deal.II/base/geometry_info.h>
30 #include <deal.II/dofs/dof_handler.h>
31 #include <deal.II/dofs/dof_accessor.h>
32 #include <deal.II/fe/fe.h>
33 #include <deal.II/fe/fe_values.h>
34 #include <deal.II/hp/fe_values.h>
35 #include <deal.II/fe/fe_update_flags.h>
36 #include <deal.II/fe/mapping_q1.h>
37 #include <deal.II/hp/q_collection.h>
38 #include <deal.II/hp/mapping_collection.h>
39 #include <deal.II/numerics/error_estimator.h>
40 #include <deal.II/distributed/tria.h>
41 
42 #include <deal.II/base/std_cxx11/bind.h>
43 
44 #include <numeric>
45 #include <algorithm>
46 #include <cmath>
47 #include <vector>
48 
49 DEAL_II_NAMESPACE_OPEN
50 
51 
52 
53 template <int spacedim>
54 template <typename InputVector, typename DoFHandlerType>
55 void
57 estimate (const Mapping<1,spacedim> &mapping,
58  const DoFHandlerType &dof_handler,
59  const Quadrature<0> &quadrature,
60  const typename FunctionMap<spacedim>::type &neumann_bc,
61  const InputVector &solution,
62  Vector<float> &error,
63  const ComponentMask &component_mask,
64  const Function<spacedim> *coefficients,
65  const unsigned int n_threads,
66  const types::subdomain_id subdomain_id,
67  const types::material_id material_id)
68 {
69  // just pass on to the other function
70  const std::vector<const InputVector *> solutions (1, &solution);
71  std::vector<Vector<float>*> errors (1, &error);
72  estimate (mapping, dof_handler, quadrature, neumann_bc, solutions, errors,
73  component_mask, coefficients, n_threads, subdomain_id, material_id);
74 }
75 
76 
77 
78 template <int spacedim>
79 template <typename InputVector, typename DoFHandlerType>
80 void
82 estimate (const DoFHandlerType &dof_handler,
83  const Quadrature<0> &quadrature,
84  const typename FunctionMap<spacedim>::type &neumann_bc,
85  const InputVector &solution,
86  Vector<float> &error,
87  const ComponentMask &component_mask,
88  const Function<spacedim> *coefficients,
89  const unsigned int n_threads,
90  const types::subdomain_id subdomain_id,
91  const types::material_id material_id)
92 {
93  estimate(StaticMappingQ1<1,spacedim>::mapping, dof_handler, quadrature, neumann_bc, solution,
94  error, component_mask, coefficients, n_threads, subdomain_id, material_id);
95 }
96 
97 
98 
99 template <int spacedim>
100 template <typename InputVector, typename DoFHandlerType>
101 void
103 estimate (const DoFHandlerType &dof_handler,
104  const Quadrature<0> &quadrature,
105  const typename FunctionMap<spacedim>::type &neumann_bc,
106  const std::vector<const InputVector *> &solutions,
107  std::vector<Vector<float>*> &errors,
108  const ComponentMask &component_mask,
109  const Function<spacedim> *coefficients,
110  const unsigned int n_threads,
111  const types::subdomain_id subdomain_id,
112  const types::material_id material_id)
113 {
114  estimate(StaticMappingQ1<1,spacedim>::mapping, dof_handler, quadrature, neumann_bc, solutions,
115  errors, component_mask, coefficients, n_threads, subdomain_id, material_id);
116 }
117 
118 
119 
120 template <int spacedim>
121 template <typename InputVector, typename DoFHandlerType>
122 void
125  const DoFHandlerType &dof_handler,
126  const hp::QCollection<0> &quadrature,
127  const typename FunctionMap<spacedim>::type &neumann_bc,
128  const InputVector &solution,
129  Vector<float> &error,
130  const ComponentMask &component_mask,
131  const Function<spacedim> *coefficients,
132  const unsigned int n_threads,
133  const types::subdomain_id subdomain_id,
134  const types::material_id material_id)
135 {
136  // just pass on to the other function
137  const std::vector<const InputVector *> solutions (1, &solution);
138  std::vector<Vector<float>*> errors (1, &error);
139  estimate (mapping, dof_handler, quadrature, neumann_bc, solutions, errors,
140  component_mask, coefficients, n_threads, subdomain_id, material_id);
141 }
142 
143 
144 template <int spacedim>
145 template <typename InputVector, typename DoFHandlerType>
146 void
148 estimate (const DoFHandlerType &dof_handler,
149  const hp::QCollection<0> &quadrature,
150  const typename FunctionMap<spacedim>::type &neumann_bc,
151  const InputVector &solution,
152  Vector<float> &error,
153  const ComponentMask &component_mask,
154  const Function<spacedim> *coefficients,
155  const unsigned int n_threads,
156  const types::subdomain_id subdomain_id,
157  const types::material_id material_id)
158 {
159  estimate(StaticMappingQ1<1,spacedim>::mapping, dof_handler, quadrature, neumann_bc, solution,
160  error, component_mask, coefficients, n_threads, subdomain_id, material_id);
161 }
162 
163 
164 
165 template <int spacedim>
166 template <typename InputVector, typename DoFHandlerType>
167 void
169 estimate (const DoFHandlerType &dof_handler,
170  const hp::QCollection<0> &quadrature,
171  const typename FunctionMap<spacedim>::type &neumann_bc,
172  const std::vector<const InputVector *> &solutions,
173  std::vector<Vector<float>*> &errors,
174  const ComponentMask &component_mask,
175  const Function<spacedim> *coefficients,
176  const unsigned int n_threads,
177  const types::subdomain_id subdomain_id,
178  const types::material_id material_id)
179 {
180  estimate(StaticMappingQ1<1,spacedim>::mapping, dof_handler, quadrature, neumann_bc, solutions,
181  errors, component_mask, coefficients, n_threads, subdomain_id, material_id);
182 }
183 
184 
185 
186 
187 template <int spacedim>
188 template <typename InputVector, typename DoFHandlerType>
190 estimate (const Mapping<1,spacedim> & /*mapping*/,
191  const DoFHandlerType & /*dof_handler*/,
192  const hp::QCollection<0> &,
193  const typename FunctionMap<spacedim>::type & /*neumann_bc*/,
194  const std::vector<const InputVector *> & /*solutions*/,
195  std::vector<Vector<float>*> & /*errors*/,
196  const ComponentMask & /*component_mask_*/,
197  const Function<spacedim> * /*coefficient*/,
198  const unsigned int,
199  const types::subdomain_id /*subdomain_id*/,
200  const types::material_id /*material_id*/)
201 {
202  Assert (false, ExcInternalError());
203 }
204 
205 
206 
207 template <int spacedim>
208 template <typename InputVector, typename DoFHandlerType>
211  const DoFHandlerType &dof_handler,
212  const Quadrature<0> &,
213  const typename FunctionMap<spacedim>::type &neumann_bc,
214  const std::vector<const InputVector *> &solutions,
215  std::vector<Vector<float>*> &errors,
216  const ComponentMask &component_mask,
217  const Function<spacedim> *coefficient,
218  const unsigned int,
219  const types::subdomain_id subdomain_id_,
220  const types::material_id material_id)
221 {
222 #ifdef DEAL_II_WITH_P4EST
223  if (dynamic_cast<const parallel::distributed::Triangulation<1,spacedim>*>
224  (&dof_handler.get_triangulation())
225  != 0)
226  Assert ((subdomain_id_ == numbers::invalid_subdomain_id)
227  ||
228  (subdomain_id_ ==
230  (dof_handler.get_triangulation()).locally_owned_subdomain()),
231  ExcMessage ("For parallel distributed triangulations, the only "
232  "valid subdomain_id that can be passed here is the "
233  "one that corresponds to the locally owned subdomain id."));
234 
235  const types::subdomain_id subdomain_id
236  = ((dynamic_cast<const parallel::distributed::Triangulation<1,spacedim>*>
237  (&dof_handler.get_triangulation())
238  != 0)
239  ?
241  (dof_handler.get_triangulation()).locally_owned_subdomain()
242  :
243  subdomain_id_);
244 #else
245  const types::subdomain_id subdomain_id
246  = subdomain_id_;
247 #endif
248 
249  const unsigned int n_components = dof_handler.get_fe().n_components();
250  const unsigned int n_solution_vectors = solutions.size();
251 
252  // sanity checks
253  Assert (neumann_bc.find(numbers::internal_face_boundary_id) == neumann_bc.end(),
254  ExcMessage("You are not allowed to list the special boundary "
255  "indicator for internal boundaries in your boundary "
256  "value map."));
257 
258  for (typename FunctionMap<spacedim>::type::const_iterator i=neumann_bc.begin();
259  i!=neumann_bc.end(); ++i)
260  Assert (i->second->n_components == n_components,
261  ExcInvalidBoundaryFunction(i->first,
262  i->second->n_components,
263  n_components));
264 
265  Assert (component_mask.represents_n_components(n_components),
266  ExcInvalidComponentMask());
267  Assert (component_mask.n_selected_components(n_components) > 0,
268  ExcInvalidComponentMask());
269 
270  Assert ((coefficient == 0) ||
271  (coefficient->n_components == n_components) ||
272  (coefficient->n_components == 1),
273  ExcInvalidCoefficient());
274 
275  Assert (solutions.size() > 0,
276  ExcNoSolutions());
277  Assert (solutions.size() == errors.size(),
278  ExcIncompatibleNumberOfElements(solutions.size(), errors.size()));
279  for (unsigned int n=0; n<solutions.size(); ++n)
280  Assert (solutions[n]->size() == dof_handler.n_dofs(),
281  ExcDimensionMismatch(solutions[n]->size(),
282  dof_handler.n_dofs()));
283 
284  Assert ((coefficient == 0) ||
285  (coefficient->n_components == n_components) ||
286  (coefficient->n_components == 1),
287  ExcInvalidCoefficient());
288 
289  for (typename FunctionMap<spacedim>::type::const_iterator i=neumann_bc.begin();
290  i!=neumann_bc.end(); ++i)
291  Assert (i->second->n_components == n_components,
292  ExcInvalidBoundaryFunction(i->first,
293  i->second->n_components,
294  n_components));
295 
296  // reserve one slot for each cell and set it to zero
297  for (unsigned int n=0; n<n_solution_vectors; ++n)
298  (*errors[n]).reinit (dof_handler.get_triangulation().n_active_cells());
299 
300  // fields to get the gradients on the present and the neighbor cell.
301  //
302  // for the neighbor gradient, we need several auxiliary fields, depending on
303  // the way we get it (see below)
304  std::vector<std::vector<std::vector<Tensor<1,spacedim,typename InputVector::value_type> > > >
305  gradients_here (n_solution_vectors,
306  std::vector<std::vector<Tensor<1,spacedim,typename InputVector::value_type> > >(2, std::vector<Tensor<1,spacedim,typename InputVector::value_type> >(n_components)));
307  std::vector<std::vector<std::vector<Tensor<1,spacedim,typename InputVector::value_type> > > >
308  gradients_neighbor (gradients_here);
309  std::vector<Vector<typename InputVector::value_type> >
310  grad_neighbor (n_solution_vectors, Vector<typename InputVector::value_type>(n_components));
311 
312  // reserve some space for coefficient values at one point. if there is no
313  // coefficient, then we fill it by unity once and for all and don't set it
314  // any more
315  Vector<double> coefficient_values (n_components);
316  if (coefficient == 0)
317  for (unsigned int c=0; c<n_components; ++c)
318  coefficient_values(c) = 1;
319 
320  const QTrapez<1> quadrature;
321  const hp::QCollection<1> q_collection(quadrature);
322  const QGauss<0> face_quadrature(1);
323  const hp::QCollection<0> q_face_collection(face_quadrature);
324 
325  const hp::FECollection<1,spacedim> fe (dof_handler.get_fe());
326 
327  hp::MappingCollection<1,spacedim> mapping_collection;
328  mapping_collection.push_back (mapping);
329 
330  hp::FEValues<1,spacedim> fe_values (mapping_collection, fe, q_collection,
332  hp::FEFaceValues<1,spacedim> fe_face_values (/*mapping_collection,*/ fe, q_face_collection,
334 
335  // loop over all cells and do something on the cells which we're told to
336  // work on. note that the error indicator is only a sum over the two
337  // contributions from the two vertices of each cell.
338  for (typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active();
339  cell != dof_handler.end();
340  ++cell)
341  if (((subdomain_id == numbers::invalid_subdomain_id)
342  ||
343  (cell->subdomain_id() == subdomain_id))
344  &&
345  ((material_id == numbers::invalid_material_id)
346  ||
347  (cell->material_id() == material_id)))
348  {
349  for (unsigned int n=0; n<n_solution_vectors; ++n)
350  (*errors[n])(cell->active_cell_index()) = 0;
351 
352  fe_values.reinit (cell);
353  for (unsigned int s=0; s<n_solution_vectors; ++s)
354  fe_values.get_present_fe_values()
355  .get_function_gradients (*solutions[s], gradients_here[s]);
356 
357  // loop over the two points bounding this line. n==0 is left point,
358  // n==1 is right point
359  for (unsigned int n=0; n<2; ++n)
360  {
361  // find left or right active neighbor
362  typename DoFHandlerType::cell_iterator neighbor = cell->neighbor(n);
363  if (neighbor.state() == IteratorState::valid)
364  while (neighbor->has_children())
365  neighbor = neighbor->child(n==0 ? 1 : 0);
366 
367  fe_face_values.reinit (cell, n);
368  Tensor<1,spacedim> normal =
369  fe_face_values.get_present_fe_values().get_all_normal_vectors()[0];
370 
371  if (neighbor.state() == IteratorState::valid)
372  {
373  fe_values.reinit (neighbor);
374 
375  for (unsigned int s=0; s<n_solution_vectors; ++s)
376  fe_values.get_present_fe_values()
377  .get_function_gradients (*solutions[s],
378  gradients_neighbor[s]);
379 
380  fe_face_values.reinit (neighbor, n==0 ? 1 : 0);
381  Tensor<1,spacedim> neighbor_normal =
382  fe_face_values.get_present_fe_values().get_all_normal_vectors()[0];
383 
384  // extract the gradient in normal direction of all the components.
385  for (unsigned int s=0; s<n_solution_vectors; ++s)
386  for (unsigned int c=0; c<n_components; ++c)
387  grad_neighbor[s](c)
388  = - (gradients_neighbor[s][n==0 ? 1 : 0][c]*neighbor_normal);
389  }
390  else if (neumann_bc.find(n) != neumann_bc.end())
391  // if Neumann b.c., then fill the gradients field which will be
392  // used later on.
393  {
394  if (n_components==1)
395  {
396  const double
397  v = neumann_bc.find(n)->second->value(cell->vertex(n));
398 
399  for (unsigned int s=0; s<n_solution_vectors; ++s)
400  grad_neighbor[s](0) = v;
401  }
402  else
403  {
404  Vector<double> v(n_components);
405  neumann_bc.find(n)->second->vector_value(cell->vertex(n), v);
406 
407  for (unsigned int s=0; s<n_solution_vectors; ++s)
408  grad_neighbor[s] = v;
409  }
410  }
411  else
412  // fill with zeroes.
413  for (unsigned int s=0; s<n_solution_vectors; ++s)
414  grad_neighbor[s] = 0;
415 
416  // if there is a coefficient, then evaluate it at the present
417  // position. if there is none, reuse the preset values.
418  if (coefficient != 0)
419  {
420  if (coefficient->n_components == 1)
421  {
422  const double c_value = coefficient->value (cell->vertex(n));
423  for (unsigned int c=0; c<n_components; ++c)
424  coefficient_values(c) = c_value;
425  }
426  else
427  coefficient->vector_value(cell->vertex(n),
428  coefficient_values);
429  }
430 
431 
432  for (unsigned int s=0; s<n_solution_vectors; ++s)
433  for (unsigned int component=0; component<n_components; ++component)
434  if (component_mask[component] == true)
435  {
436  // get gradient here
437  const double grad_here = gradients_here[s][n][component]
438  * normal;
439 
440  const double jump = ((grad_here - grad_neighbor[s](component)) *
441  coefficient_values(component));
442  (*errors[s])(cell->active_cell_index()) += jump*jump * cell->diameter();
443  }
444  }
445 
446  for (unsigned int s=0; s<n_solution_vectors; ++s)
447  (*errors[s])(cell->active_cell_index()) = std::sqrt((*errors[s])(cell->active_cell_index()));
448  }
449 }
450 
451 
452 // explicit instantiations
453 #include "error_estimator_1d.inst"
454 
455 
456 DEAL_II_NAMESPACE_CLOSE
const types::subdomain_id invalid_subdomain_id
Definition: types.h:239
unsigned char material_id
Definition: types.h:130
::ExceptionBase & ExcMessage(std::string arg1)
bool represents_n_components(const unsigned int n) const
const ::FEValues< dim, spacedim > & get_present_fe_values() const
Definition: fe_values.h:598
const unsigned int n_components
Definition: function.h:144
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, 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:144
void push_back(const FiniteElement< dim, spacedim > &new_fe)
#define Assert(cond, exc)
Definition: exceptions.h:294
Abstract base class for mapping classes.
Definition: dof_tools.h:52
std::map< types::boundary_id, const Function< dim, Number > * > type
Definition: function_map.h:81
unsigned int subdomain_id
Definition: types.h:42
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, 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:261
virtual Number value(const Point< dim > &p, const unsigned int component=0) const
Definition: mpi.h:48
Shape function gradients.
Normal vectors.
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
Iterator points to a valid object.
const types::boundary_id internal_face_boundary_id
Definition: types.h:210
const types::material_id invalid_material_id
Definition: types.h:185
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< dim-1 > &quadrature, const typename FunctionMap< spacedim >::type &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=0, 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)
virtual void vector_value(const Point< dim > &p, Vector< Number > &values) const