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> 43 #include <deal.II/base/std_cxx11/bind.h> 50 DEAL_II_NAMESPACE_OPEN
54 template <
int spacedim>
55 template <
typename InputVector,
typename DoFHandlerType>
59 const DoFHandlerType &dof_handler,
62 const InputVector &solution,
66 const unsigned int n_threads,
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);
79 template <
int spacedim>
80 template <
typename InputVector,
typename DoFHandlerType>
86 const InputVector &solution,
90 const unsigned int n_threads,
95 error, component_mask, coefficients, n_threads, subdomain_id, material_id);
100 template <
int spacedim>
101 template <
typename InputVector,
typename DoFHandlerType>
107 const std::vector<const InputVector *> &solutions,
111 const unsigned int n_threads,
116 errors, component_mask, coefficients, n_threads, subdomain_id, material_id);
121 template <
int spacedim>
122 template <
typename InputVector,
typename DoFHandlerType>
126 const DoFHandlerType &dof_handler,
129 const InputVector &solution,
133 const unsigned int n_threads,
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);
145 template <
int spacedim>
146 template <
typename InputVector,
typename DoFHandlerType>
152 const InputVector &solution,
156 const unsigned int n_threads,
161 error, component_mask, coefficients, n_threads, subdomain_id, material_id);
166 template <
int spacedim>
167 template <
typename InputVector,
typename DoFHandlerType>
173 const std::vector<const InputVector *> &solutions,
177 const unsigned int n_threads,
182 errors, component_mask, coefficients, n_threads, subdomain_id, material_id);
188 template <
int spacedim>
189 template <
typename InputVector,
typename DoFHandlerType>
192 const DoFHandlerType & ,
195 const std::vector<const InputVector *> & ,
208 template <
int spacedim>
209 template <
typename InputVector,
typename DoFHandlerType>
212 const DoFHandlerType &dof_handler,
215 const std::vector<const InputVector *> &solutions,
223 typedef typename InputVector::value_type number;
224 #ifdef DEAL_II_WITH_P4EST 226 (&dof_handler.get_triangulation())
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."));
239 (&dof_handler.get_triangulation())
243 (dof_handler.get_triangulation()).locally_owned_subdomain()
251 const unsigned int n_components = dof_handler.get_fe().n_components();
252 const unsigned int n_solution_vectors = solutions.size();
256 ExcMessage(
"You are not allowed to list the special boundary " 257 "indicator for internal boundaries in your boundary " 261 i!=neumann_bc.end(); ++i)
262 Assert (i->second->n_components == n_components,
264 i->second->n_components,
272 Assert ((coefficient == 0) ||
277 Assert (solutions.size() > 0,
279 Assert (solutions.size() == errors.size(),
281 for (
unsigned int n=0; n<solutions.size(); ++n)
282 Assert (solutions[n]->size() == dof_handler.n_dofs(),
284 dof_handler.n_dofs()));
286 Assert ((coefficient == 0) ||
292 i!=neumann_bc.end(); ++i)
293 Assert (i->second->n_components == n_components,
295 i->second->n_components,
299 for (
unsigned int n=0; n<n_solution_vectors; ++n)
300 (*errors[n]).reinit (dof_handler.get_triangulation().n_active_cells());
306 std::vector<std::vector<std::vector<Tensor<1,spacedim,number> > > >
307 gradients_here (n_solution_vectors,
309 std::vector<std::vector<std::vector<Tensor<1,spacedim,number> > > >
310 gradients_neighbor (gradients_here);
311 std::vector<Vector<number> >
318 if (coefficient == 0)
319 for (
unsigned int c=0; c<n_components; ++c)
320 coefficient_values(c) = 1;
340 for (
typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active();
341 cell != dof_handler.end();
345 (cell->subdomain_id() == subdomain_id))
349 (cell->material_id() == material_id)))
351 for (
unsigned int n=0; n<n_solution_vectors; ++n)
352 (*errors[n])(cell->active_cell_index()) = 0;
355 for (
unsigned int s=0; s<n_solution_vectors; ++s)
357 .get_function_gradients (*solutions[s], gradients_here[s]);
361 for (
unsigned int n=0; n<2; ++n)
364 typename DoFHandlerType::cell_iterator neighbor = cell->neighbor(n);
366 while (neighbor->has_children())
367 neighbor = neighbor->child(n==0 ? 1 : 0);
369 fe_face_values.
reinit (cell, n);
375 fe_values.
reinit (neighbor);
377 for (
unsigned int s=0; s<n_solution_vectors; ++s)
379 .get_function_gradients (*solutions[s],
380 gradients_neighbor[s]);
382 fe_face_values.
reinit (neighbor, n==0 ? 1 : 0);
387 for (
unsigned int s=0; s<n_solution_vectors; ++s)
388 for (
unsigned int c=0; c<n_components; ++c)
390 = - (gradients_neighbor[s][n==0 ? 1 : 0][c]*neighbor_normal);
392 else if (neumann_bc.find(n) != neumann_bc.end())
399 v = neumann_bc.find(n)->second->value(cell->vertex(n));
401 for (
unsigned int s=0; s<n_solution_vectors; ++s)
402 grad_neighbor[s](0) = v;
407 neumann_bc.find(n)->second->vector_value(cell->vertex(n), v);
409 for (
unsigned int s=0; s<n_solution_vectors; ++s)
410 grad_neighbor[s] = v;
415 for (
unsigned int s=0; s<n_solution_vectors; ++s)
416 grad_neighbor[s] = 0;
420 if (coefficient != 0)
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;
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)
439 const number grad_here = gradients_here[s][n][component]
442 const number jump = ((grad_here - grad_neighbor[s](component)) *
443 coefficient_values(component));
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()));
455 #include "error_estimator_1d.inst" 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
static ::ExceptionBase & ExcInvalidComponentMask()
unsigned char material_id
bool represents_n_components(const unsigned int n) const
const ::FEValues< dim, spacedim > & get_present_fe_values() const
const unsigned int n_components
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)
static ::ExceptionBase & ExcMessage(std::string arg1)
void push_back(const FiniteElement< dim, spacedim > &new_fe)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
std::map< types::boundary_id, const Function< dim, Number > * > type
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
unsigned int subdomain_id
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)
virtual Number value(const Point< dim > &p, const unsigned int component=0) const
Shape function gradients.
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
const types::material_id invalid_material_id
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()