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