16 #include <deal.II/base/point.h> 17 #include <deal.II/base/auto_derivative_function.h> 18 #include <deal.II/lac/vector.h> 22 DEAL_II_NAMESPACE_OPEN
28 const double initial_time)
30 Function<dim>(n_components, initial_time),
59 for (
unsigned int i=0; i<dim; ++i)
67 const unsigned int comp)
const 75 for (
unsigned int i=0; i<dim; ++i)
78 grad[i]=(this->
value(p, comp)-this->
value(q1, comp))/
h;
85 for (
unsigned int i=0; i<dim; ++i)
89 grad[i]=(this->
value(q1, comp)-this->
value(q2, comp))/(2*
h);
96 for (
unsigned int i=0; i<dim; ++i)
102 grad[i]=(- this->
value(q1, comp)
103 +8*this->
value(q2, comp)
104 -8*this->
value(q3, comp)
105 + this->
value(q4, comp))/(12*
h);
111 Assert(
false, ExcInvalidFormula());
124 ExcDimensionMismatch(gradients.size(), this->
n_components));
132 const double h_inv=1./
h;
133 for (
unsigned int i=0; i<dim; ++i)
139 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
140 gradients[comp][i]=(v(comp)-v1(comp))*h_inv;
149 const double h_inv_2=1./(2*
h);
150 for (
unsigned int i=0; i<dim; ++i)
157 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
158 gradients[comp][i]=(v1(comp)-v2(comp))*h_inv_2;
169 const double h_inv_12=1./(12*
h);
170 for (
unsigned int i=0; i<dim; ++i)
181 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
182 gradients[comp][i]=(-v1(comp)+8*v2(comp)-8*v3(comp)+v4(comp))*h_inv_12;
188 Assert(
false, ExcInvalidFormula());
198 const unsigned int comp)
const 200 Assert (gradients.size() == points.size(),
201 ExcDimensionMismatch(gradients.size(), points.size()));
208 for (
unsigned int p=0; p<points.size(); ++p)
209 for (
unsigned int i=0; i<dim; ++i)
212 gradients[p][i]=(this->
value(points[p], comp)-this->
value(q1, comp))/
h;
220 for (
unsigned int p=0; p<points.size(); ++p)
221 for (
unsigned int i=0; i<dim; ++i)
225 gradients[p][i]=(this->
value(q1, comp)-this->
value(q2, comp))/(2*
h);
233 for (
unsigned int p=0; p<points.size(); ++p)
234 for (
unsigned int i=0; i<dim; ++i)
240 gradients[p][i]=(- this->
value(q1, comp)
241 +8*this->
value(q2, comp)
242 -8*this->
value(q3, comp)
243 + this->
value(q4, comp))/(12*
h);
249 Assert(
false, ExcInvalidFormula());
261 Assert (gradients.size() == points.size(),
262 ExcDimensionMismatch(gradients.size(), points.size()));
263 for (
unsigned int p=0; p<points.size(); ++p)
265 ExcDimensionMismatch(gradients.size(), this->
n_components));
272 for (
unsigned int p=0; p<points.size(); ++p)
273 for (
unsigned int i=0; i<dim; ++i)
276 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
277 gradients[p][comp][i]=(this->
value(points[p], comp)-this->
value(q1, comp))/
h;
285 for (
unsigned int p=0; p<points.size(); ++p)
286 for (
unsigned int i=0; i<dim; ++i)
290 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
291 gradients[p][comp][i]=(this->
value(q1, comp) -
292 this->
value(q2, comp))/(2*
h);
300 for (
unsigned int p=0; p<points.size(); ++p)
301 for (
unsigned int i=0; i<dim; ++i)
307 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
308 gradients[p][comp][i]=(- this->
value(q1, comp)
309 +8*this->
value(q2, comp)
310 -8*this->
value(q3, comp)
311 + this->
value(q4, comp))/(12*
h);
317 Assert(
false, ExcInvalidFormula());
337 Assert(
false, ExcNotImplemented());
347 DEAL_II_NAMESPACE_CLOSE
virtual void vector_gradient(const Point< dim > &p, std::vector< Tensor< 1, dim > > &gradients) const
static DifferenceFormula get_formula_of_order(const unsigned int ord)
DifferenceFormula formula
virtual void gradient_list(const std::vector< Point< dim > > &points, std::vector< Tensor< 1, dim > > &gradients, const unsigned int component=0) const
virtual ~AutoDerivativeFunction()
void set_formula(const DifferenceFormula formula=Euler)
const unsigned int n_components
#define Assert(cond, exc)
void set_h(const double h)
virtual void vector_gradient_list(const std::vector< Point< dim > > &points, std::vector< std::vector< Tensor< 1, dim > > > &gradients) const
AutoDerivativeFunction(const double h, const unsigned int n_components=1, const double initial_time=0.0)
virtual double value(const Point< dim > &p, const unsigned int component=0) const
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
std::vector< Tensor< 1, dim > > ht
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual void vector_value(const Point< dim > &p, Vector< double > &values) const