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)
60 "match any known difference formula."));
72 for (
unsigned int i=0; i<dim; ++i)
80 const unsigned int comp)
const 88 for (
unsigned int i=0; i<dim; ++i)
91 grad[i]=(this->value(p, comp)-this->value(q1, comp))/h;
98 for (
unsigned int i=0; i<dim; ++i)
102 grad[i]=(this->value(q1, comp)-this->value(q2, comp))/(2*h);
109 for (
unsigned int i=0; i<dim; ++i)
115 grad[i]=(- this->value(q1, comp)
116 +8*this->value(q2, comp)
117 -8*this->value(q3, comp)
118 + this->value(q4, comp))/(12*h);
145 const double h_inv=1./h;
146 for (
unsigned int i=0; i<dim; ++i)
149 this->vector_value(p, v);
150 this->vector_value(q1, v1);
152 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
153 gradients[comp][i]=(v(comp)-v1(comp))*h_inv;
162 const double h_inv_2=1./(2*h);
163 for (
unsigned int i=0; i<dim; ++i)
167 this->vector_value(q1, v1);
168 this->vector_value(q2, v2);
170 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
171 gradients[comp][i]=(v1(comp)-v2(comp))*h_inv_2;
182 const double h_inv_12=1./(12*h);
183 for (
unsigned int i=0; i<dim; ++i)
189 this->vector_value(q1, v1);
190 this->vector_value(q2, v2);
191 this->vector_value(q3, v3);
192 this->vector_value(q4, v4);
194 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
195 gradients[comp][i]=(-v1(comp)+8*v2(comp)-8*v3(comp)+v4(comp))*h_inv_12;
211 const unsigned int comp)
const 213 Assert (gradients.size() == points.size(),
221 for (
unsigned int p=0; p<points.size(); ++p)
222 for (
unsigned int i=0; i<dim; ++i)
225 gradients[p][i]=(this->value(points[p], comp)-this->value(q1, comp))/h;
233 for (
unsigned int p=0; p<points.size(); ++p)
234 for (
unsigned int i=0; i<dim; ++i)
238 gradients[p][i]=(this->value(q1, comp)-this->value(q2, comp))/(2*h);
246 for (
unsigned int p=0; p<points.size(); ++p)
247 for (
unsigned int i=0; i<dim; ++i)
253 gradients[p][i]=(- this->value(q1, comp)
254 +8*this->value(q2, comp)
255 -8*this->value(q3, comp)
256 + this->value(q4, comp))/(12*h);
274 Assert (gradients.size() == points.size(),
276 for (
unsigned int p=0; p<points.size(); ++p)
285 for (
unsigned int p=0; p<points.size(); ++p)
286 for (
unsigned int i=0; i<dim; ++i)
289 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
290 gradients[p][comp][i]=(this->value(points[p], comp)-this->value(q1, comp))/h;
298 for (
unsigned int p=0; p<points.size(); ++p)
299 for (
unsigned int i=0; i<dim; ++i)
303 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
304 gradients[p][comp][i]=(this->value(q1, comp) -
305 this->value(q2, comp))/(2*h);
313 for (
unsigned int p=0; p<points.size(); ++p)
314 for (
unsigned int i=0; i<dim; ++i)
320 for (
unsigned int comp=0; comp<this->
n_components; ++comp)
321 gradients[p][comp][i]=(- this->value(q1, comp)
322 +8*this->value(q2, comp)
323 -8*this->value(q3, comp)
324 + this->value(q4, comp))/(12*h);
360 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)
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)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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 Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
static ::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)