Reference documentation for deal.II version 8.5.1
fe_field_function.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2007 - 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 #ifndef dealii__fe_function_h
17 #define dealii__fe_function_h
18 
19 #include <deal.II/base/function.h>
20 #include <deal.II/dofs/dof_handler.h>
21 #include <deal.II/dofs/dof_accessor.h>
22 #include <deal.II/fe/mapping_q1.h>
23 #include <deal.II/base/function.h>
24 #include <deal.II/base/point.h>
25 #include <deal.II/base/tensor.h>
26 #include <deal.II/base/thread_local_storage.h>
27 
28 #include <deal.II/lac/vector.h>
29 
30 #include <boost/optional.hpp>
31 
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 namespace VectorTools
36 {
38 }
39 
40 namespace Functions
41 {
42 
162  template <int dim,
163  typename DoFHandlerType=DoFHandler<dim>,
164  typename VectorType=Vector<double> >
165  class FEFieldFunction : public Function<dim, typename VectorType::value_type>
166  {
167  public:
176  FEFieldFunction (const DoFHandlerType &dh,
177  const VectorType &data_vector,
179 
185  void set_active_cell (const typename DoFHandlerType::active_cell_iterator &newcell);
186 
202  virtual void vector_value (const Point<dim> &p,
204 
222  virtual typename VectorType::value_type value (const Point< dim > &p,
223  const unsigned int component = 0) const;
224 
240  virtual void value_list (const std::vector<Point< dim > > &points,
241  std::vector<typename VectorType::value_type > &values,
242  const unsigned int component = 0) const;
243 
244 
260  virtual void vector_value_list (const std::vector<Point< dim > > &points,
261  std::vector<Vector<typename VectorType::value_type> > &values) const;
262 
278  virtual void
279  vector_gradient (const Point< dim > &p,
280  std::vector< Tensor< 1, dim,typename VectorType::value_type > > &gradients) const;
281 
298  const unsigned int component = 0)const;
299 
313  virtual void
314  vector_gradient_list (const std::vector< Point< dim > > &p,
315  std::vector<
316  std::vector< Tensor< 1, dim,typename VectorType::value_type > > > &gradients) const;
317 
331  virtual void
332  gradient_list (const std::vector< Point< dim > > &p,
333  std::vector< Tensor< 1, dim,typename VectorType::value_type > > &gradients,
334  const unsigned int component=0) const;
335 
336 
348  virtual typename VectorType::value_type
349  laplacian (const Point<dim> &p,
350  const unsigned int component = 0) const;
351 
364  virtual void
365  vector_laplacian (const Point<dim> &p,
367 
379  virtual void
380  laplacian_list (const std::vector<Point<dim> > &points,
381  std::vector<typename VectorType::value_type> &values,
382  const unsigned int component = 0) const;
383 
395  virtual void
396  vector_laplacian_list (const std::vector<Point<dim> > &points,
397  std::vector<Vector<typename VectorType::value_type> > &values) const;
398 
410  unsigned int
412  (const std::vector<Point<dim> > &points,
413  std::vector<typename DoFHandlerType::active_cell_iterator > &cells,
414  std::vector<std::vector<Point<dim> > > &qpoints,
415  std::vector<std::vector<unsigned int> > &maps) const;
416 
421 
422  private:
426  typedef
429 
434 
438  const VectorType &data_vector;
439 
444 
449 
455  boost::optional<Point<dim> >
456  get_reference_coordinates (const typename DoFHandlerType::active_cell_iterator &cell,
457  const Point<dim> &point) const;
458  };
459 }
460 
461 DEAL_II_NAMESPACE_CLOSE
462 
463 #endif
SmartPointer< const DoFHandlerType, FEFieldFunction< dim, DoFHandlerType, VectorType > > dh
FEFieldFunction(const DoFHandlerType &dh, const VectorType &data_vector, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
virtual VectorType::value_type value(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value(const Point< dim > &p, Vector< typename VectorType::value_type > &values) const
virtual void gradient_list(const std::vector< Point< dim > > &p, std::vector< Tensor< 1, dim, typename VectorType::value_type > > &gradients, const unsigned int component=0) const
virtual void vector_laplacian(const Point< dim > &p, Vector< typename VectorType::value_type > &values) const
virtual void laplacian_list(const std::vector< Point< dim > > &points, std::vector< typename VectorType::value_type > &values, const unsigned int component=0) const
const Mapping< dim > & mapping
void set_active_cell(const typename DoFHandlerType::active_cell_iterator &newcell)
boost::optional< Point< dim > > get_reference_coordinates(const typename DoFHandlerType::active_cell_iterator &cell, const Point< dim > &point) const
virtual VectorType::value_type laplacian(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value_list(const std::vector< Point< dim > > &points, std::vector< Vector< typename VectorType::value_type > > &values) const
static ::ExceptionBase & ExcPointNotAvailableHere()
virtual void vector_gradient_list(const std::vector< Point< dim > > &p, std::vector< std::vector< Tensor< 1, dim, typename VectorType::value_type > > > &gradients) const
VectorTools::ExcPointNotAvailableHere ExcPointNotAvailableHere
const VectorType & data_vector
virtual void vector_laplacian_list(const std::vector< Point< dim > > &points, std::vector< Vector< typename VectorType::value_type > > &values) const
Definition: mpi.h:41
Threads::ThreadLocalStorage< typename DoFHandlerType::active_cell_iterator > cell_hint_t
virtual void vector_gradient(const Point< dim > &p, std::vector< Tensor< 1, dim, typename VectorType::value_type > > &gradients) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< typename VectorType::value_type > &values, const unsigned int component=0) const
unsigned int compute_point_locations(const std::vector< Point< dim > > &points, std::vector< typename DoFHandlerType::active_cell_iterator > &cells, std::vector< std::vector< Point< dim > > > &qpoints, std::vector< std::vector< unsigned int > > &maps) const
virtual Tensor< 1, dim, typename VectorType::value_type > gradient(const Point< dim > &p, const unsigned int component=0) const