Reference documentation for deal.II version 9.0.0
fe_field_function.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2007 - 2018 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 #include <deal.II/grid/grid_tools_cache.h>
28 
29 #include <deal.II/lac/vector.h>
30 
31 #include <boost/optional.hpp>
32 
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 namespace VectorTools
37 {
39 }
40 
41 namespace Functions
42 {
43 
157  template <int dim,
158  typename DoFHandlerType=DoFHandler<dim>,
159  typename VectorType=Vector<double> >
160  class FEFieldFunction : public Function<dim, typename VectorType::value_type>
161  {
162  public:
171  FEFieldFunction (const DoFHandlerType &dh,
172  const VectorType &data_vector,
174 
180  void set_active_cell (const typename DoFHandlerType::active_cell_iterator &newcell);
181 
197  virtual void vector_value (const Point<dim> &p,
199 
217  virtual typename VectorType::value_type value (const Point< dim > &p,
218  const unsigned int component = 0) const;
219 
235  virtual void value_list (const std::vector<Point< dim > > &points,
236  std::vector<typename VectorType::value_type > &values,
237  const unsigned int component = 0) const;
238 
239 
255  virtual void vector_value_list (const std::vector<Point< dim > > &points,
256  std::vector<Vector<typename VectorType::value_type> > &values) const;
257 
273  virtual void
274  vector_gradient (const Point< dim > &p,
275  std::vector< Tensor< 1, dim,typename VectorType::value_type > > &gradients) const;
276 
293  const unsigned int component = 0)const;
294 
308  virtual void
309  vector_gradient_list (const std::vector< Point< dim > > &p,
310  std::vector<
311  std::vector< Tensor< 1, dim,typename VectorType::value_type > > > &gradients) const;
312 
326  virtual void
327  gradient_list (const std::vector< Point< dim > > &p,
328  std::vector< Tensor< 1, dim,typename VectorType::value_type > > &gradients,
329  const unsigned int component=0) const;
330 
331 
343  virtual typename VectorType::value_type
344  laplacian (const Point<dim> &p,
345  const unsigned int component = 0) const;
346 
359  virtual void
360  vector_laplacian (const Point<dim> &p,
362 
374  virtual void
375  laplacian_list (const std::vector<Point<dim> > &points,
376  std::vector<typename VectorType::value_type> &values,
377  const unsigned int component = 0) const;
378 
390  virtual void
391  vector_laplacian_list (const std::vector<Point<dim> > &points,
392  std::vector<Vector<typename VectorType::value_type> > &values) const;
393 
416  unsigned int
418  (const std::vector<Point<dim> > &points,
419  std::vector<typename DoFHandlerType::active_cell_iterator > &cells,
420  std::vector<std::vector<Point<dim> > > &qpoints,
421  std::vector<std::vector<unsigned int> > &maps) const;
422 
423  private:
427  typedef
430 
435 
439  const VectorType &data_vector;
440 
445 
450 
455 
461  boost::optional<Point<dim> >
462  get_reference_coordinates (const typename DoFHandlerType::active_cell_iterator &cell,
463  const Point<dim> &point) const;
464  };
465 }
466 
467 DEAL_II_NAMESPACE_CLOSE
468 
469 #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)
GridTools::Cache< dim, DoFHandlerType::space_dimension > cache
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
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:53
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