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> 42 #include <deal.II/base/std_cxx11/bind.h> 49 DEAL_II_NAMESPACE_OPEN
53 template <
int spacedim>
54 template <
typename InputVector,
typename DoFHandlerType>
58 const DoFHandlerType &dof_handler,
61 const InputVector &solution,
65 const unsigned int n_threads,
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);
78 template <
int spacedim>
79 template <
typename InputVector,
typename DoFHandlerType>
85 const InputVector &solution,
89 const unsigned int n_threads,
94 error, component_mask, coefficients, n_threads, subdomain_id, material_id);
99 template <
int spacedim>
100 template <
typename InputVector,
typename DoFHandlerType>
106 const std::vector<const InputVector *> &solutions,
110 const unsigned int n_threads,
115 errors, component_mask, coefficients, n_threads, subdomain_id, material_id);
120 template <
int spacedim>
121 template <
typename InputVector,
typename DoFHandlerType>
125 const DoFHandlerType &dof_handler,
128 const InputVector &solution,
132 const unsigned int n_threads,
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);
144 template <
int spacedim>
145 template <
typename InputVector,
typename DoFHandlerType>
151 const InputVector &solution,
155 const unsigned int n_threads,
160 error, component_mask, coefficients, n_threads, subdomain_id, material_id);
165 template <
int spacedim>
166 template <
typename InputVector,
typename DoFHandlerType>
172 const std::vector<const InputVector *> &solutions,
176 const unsigned int n_threads,
181 errors, component_mask, coefficients, n_threads, subdomain_id, material_id);
187 template <
int spacedim>
188 template <
typename InputVector,
typename DoFHandlerType>
191 const DoFHandlerType & ,
194 const std::vector<const InputVector *> & ,
202 Assert (
false, ExcInternalError());
207 template <
int spacedim>
208 template <
typename InputVector,
typename DoFHandlerType>
211 const DoFHandlerType &dof_handler,
214 const std::vector<const InputVector *> &solutions,
222 #ifdef DEAL_II_WITH_P4EST 224 (&dof_handler.get_triangulation())
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."));
237 (&dof_handler.get_triangulation())
241 (dof_handler.get_triangulation()).locally_owned_subdomain()
249 const unsigned int n_components = dof_handler.get_fe().n_components();
250 const unsigned int n_solution_vectors = solutions.size();
254 ExcMessage(
"You are not allowed to list the special boundary " 255 "indicator for internal boundaries in your boundary " 259 i!=neumann_bc.end(); ++i)
260 Assert (i->second->n_components == n_components,
261 ExcInvalidBoundaryFunction(i->first,
262 i->second->n_components,
266 ExcInvalidComponentMask());
268 ExcInvalidComponentMask());
270 Assert ((coefficient == 0) ||
273 ExcInvalidCoefficient());
275 Assert (solutions.size() > 0,
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()));
284 Assert ((coefficient == 0) ||
287 ExcInvalidCoefficient());
290 i!=neumann_bc.end(); ++i)
291 Assert (i->second->n_components == n_components,
292 ExcInvalidBoundaryFunction(i->first,
293 i->second->n_components,
297 for (
unsigned int n=0; n<n_solution_vectors; ++n)
298 (*errors[n]).reinit (dof_handler.get_triangulation().n_active_cells());
304 std::vector<std::vector<std::vector<Tensor<1,spacedim,typename InputVector::value_type> > > >
305 gradients_here (n_solution_vectors,
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> >
316 if (coefficient == 0)
317 for (
unsigned int c=0; c<n_components; ++c)
318 coefficient_values(c) = 1;
338 for (
typename DoFHandlerType::active_cell_iterator cell = dof_handler.begin_active();
339 cell != dof_handler.end();
343 (cell->subdomain_id() == subdomain_id))
347 (cell->material_id() == material_id)))
349 for (
unsigned int n=0; n<n_solution_vectors; ++n)
350 (*errors[n])(cell->active_cell_index()) = 0;
353 for (
unsigned int s=0; s<n_solution_vectors; ++s)
355 .get_function_gradients (*solutions[s], gradients_here[s]);
359 for (
unsigned int n=0; n<2; ++n)
362 typename DoFHandlerType::cell_iterator neighbor = cell->neighbor(n);
364 while (neighbor->has_children())
365 neighbor = neighbor->child(n==0 ? 1 : 0);
367 fe_face_values.
reinit (cell, n);
373 fe_values.
reinit (neighbor);
375 for (
unsigned int s=0; s<n_solution_vectors; ++s)
377 .get_function_gradients (*solutions[s],
378 gradients_neighbor[s]);
380 fe_face_values.
reinit (neighbor, n==0 ? 1 : 0);
385 for (
unsigned int s=0; s<n_solution_vectors; ++s)
386 for (
unsigned int c=0; c<n_components; ++c)
388 = - (gradients_neighbor[s][n==0 ? 1 : 0][c]*neighbor_normal);
390 else if (neumann_bc.find(n) != neumann_bc.end())
397 v = neumann_bc.find(n)->second->value(cell->vertex(n));
399 for (
unsigned int s=0; s<n_solution_vectors; ++s)
400 grad_neighbor[s](0) = v;
405 neumann_bc.find(n)->second->vector_value(cell->vertex(n), v);
407 for (
unsigned int s=0; s<n_solution_vectors; ++s)
408 grad_neighbor[s] = v;
413 for (
unsigned int s=0; s<n_solution_vectors; ++s)
414 grad_neighbor[s] = 0;
418 if (coefficient != 0)
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;
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)
437 const double grad_here = gradients_here[s][n][component]
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();
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()));
453 #include "error_estimator_1d.inst" 456 DEAL_II_NAMESPACE_CLOSE
const types::subdomain_id invalid_subdomain_id
unsigned char material_id
::ExceptionBase & ExcMessage(std::string arg1)
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)
void push_back(const FiniteElement< dim, spacedim > &new_fe)
#define Assert(cond, exc)
Abstract base class for mapping classes.
std::map< types::boundary_id, const Function< dim, Number > * > type
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)
virtual void vector_value(const Point< dim > &p, Vector< Number > &values) const