DOLFIN
DOLFIN C++ interface
GenericFunction.h
1 // Copyright (C) 2009 Anders Logg
2 //
3 // This file is part of DOLFIN.
4 //
5 // DOLFIN is free software: you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published by
7 // the Free Software Foundation, either version 3 of the License, or
8 // (at your option) any later version.
9 //
10 // DOLFIN is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17 //
18 // Modified by Garth N. Wells, 2009.
19 //
20 // First added: 2009-09-28
21 // Last changed: 2011-01-19
22 
23 #ifndef __GENERIC_FUNCTION_H
24 #define __GENERIC_FUNCTION_H
25 
26 #include <memory>
27 #include <Eigen/Dense>
28 #include <ufc.h>
29 #include <dolfin/common/Array.h>
30 #include <dolfin/common/Variable.h>
31 
32 namespace dolfin
33 {
34 
35  class Mesh;
36  class Cell;
37  class Point;
38  class FiniteElement;
39  class FunctionSpace;
40 
52 
53  class GenericFunction : public ufc::function, public Variable
54  {
55  public:
56 
59 
61  virtual ~GenericFunction();
62 
63  //--- Functions that must be implemented by sub-classes ---
64 
66  virtual std::size_t value_rank() const = 0;
67 
69  virtual std::size_t value_dimension(std::size_t i) const = 0;
70 
72  virtual std::vector<std::size_t> value_shape() const = 0;
73 
75  virtual void eval(Array<double>& values, const Array<double>& x,
76  const ufc::cell& cell) const;
77 
79  virtual void eval(Array<double>& values, const Array<double>& x) const;
80 
82  virtual void eval(Eigen::Ref<Eigen::VectorXd> values,
83  Eigen::Ref<const Eigen::VectorXd> x,
84  const ufc::cell& cell) const;
85 
87  virtual void eval(Eigen::Ref<Eigen::VectorXd> values,
88  Eigen::Ref<const Eigen::VectorXd> x) const;
89 
91  virtual void restrict(double* w,
92  const FiniteElement& element,
93  const Cell& dolfin_cell,
94  const double* coordinate_dofs,
95  const ufc::cell& ufc_cell) const = 0;
96 
98  virtual void compute_vertex_values(std::vector<double>& vertex_values,
99  const Mesh& mesh) const = 0;
100 
101  //--- Optional functions to be implemented by sub-classes ---
102 
104  virtual void update() const {}
105 
106  //--- Convenience functions ---
107 
109  double operator() (double x) const;
110 
112  double operator() (double x, double y) const;
113 
115  double operator() (double x, double y, double z) const;
116 
118  double operator() (const Point& p) const;
119 
121  void operator() (Array<double>& values, double x) const;
122 
124  void operator() (Array<double>& values, double x, double y) const;
125 
127  void operator() (Array<double>& values, double x, double y, double z) const;
128 
130  void operator() (Array<double>& values, const Point& p) const;
131 
133 
135  std::size_t value_size() const;
136 
137  //--- Implementation of ufc::function interface ---
138 
143  virtual void evaluate(double* values,
144  const double* coordinates,
145  const ufc::cell& cell) const override;
146 
148  virtual std::shared_ptr<const FunctionSpace> function_space() const = 0;
149 
150  protected:
151 
153  void restrict_as_ufc_function(double* w,
154  const FiniteElement& element,
155  const Cell& dolfin_cell,
156  const double* coordinate_dofs,
157  const ufc::cell& ufc_cell) const;
158 
159  };
160 
161 }
162 
163 #endif
virtual void restrict(double *w, const FiniteElement &element, const Cell &dolfin_cell, const double *coordinate_dofs, const ufc::cell &ufc_cell) const =0
Restrict function to local cell (compute expansion coefficients w)
std::size_t value_size() const
Evaluation at given point.
Definition: GenericFunction.cpp:181
Common base class for DOLFIN variables.
Definition: Variable.h:35
virtual ~GenericFunction()
Destructor.
Definition: GenericFunction.cpp:36
virtual void evaluate(double *values, const double *coordinates, const ufc::cell &cell) const override
Definition: GenericFunction.cpp:189
virtual void update() const
Update off-process ghost coefficients.
Definition: GenericFunction.h:104
Definition: adapt.h:29
Definition: Point.h:40
Definition: Array.h:41
virtual void compute_vertex_values(std::vector< double > &vertex_values, const Mesh &mesh) const =0
Compute values at all mesh vertices.
void restrict_as_ufc_function(double *w, const FiniteElement &element, const Cell &dolfin_cell, const double *coordinate_dofs, const ufc::cell &ufc_cell) const
Restrict as UFC function (by calling eval)
Definition: GenericFunction.cpp:205
virtual std::vector< std::size_t > value_shape() const =0
Return value shape.
A Cell is a MeshEntity of topological codimension 0.
Definition: Cell.h:42
double operator()(double x) const
Evaluation at given point (scalar function)
Definition: GenericFunction.cpp:71
GenericFunction()
Constructor.
Definition: GenericFunction.cpp:31
virtual std::size_t value_rank() const =0
Return value rank.
virtual std::size_t value_dimension(std::size_t i) const =0
Return value dimension for given axis.
Definition: GenericFunction.h:53
This is a wrapper for a UFC finite element (ufc::finite_element).
Definition: FiniteElement.h:35
virtual void eval(Array< double > &values, const Array< double > &x, const ufc::cell &cell) const
Evaluate at given point in given cell (deprecated)
Definition: GenericFunction.cpp:41
Definition: Mesh.h:82
virtual std::shared_ptr< const FunctionSpace > function_space() const =0
Pointer to FunctionSpace, if appropriate, otherwise NULL.