]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move get_data(), get_face_data(), get_subface_data() definitions directly to places...
authorLukas Korous <lukas.korous@gmail.com>
Wed, 12 Aug 2015 20:02:11 +0000 (22:02 +0200)
committerLukas Korous <lukas.korous@gmail.com>
Thu, 13 Aug 2015 19:06:35 +0000 (21:06 +0200)
include/deal.II/fe/fe_face.h
include/deal.II/fe/fe_poly.h
include/deal.II/fe/fe_poly.templates.h
include/deal.II/fe/fe_poly_face.h
include/deal.II/fe/fe_poly_face.templates.h
include/deal.II/fe/fe_poly_tensor.h
source/fe/fe_face.cc
source/fe/fe_poly_tensor.cc

index a4a54fcab43f3d7a550f9b2d491d18e7f166ed39..dd51b9343f52e260c28e9f959a1a6c1b7271f806 100644 (file)
@@ -231,17 +231,47 @@ protected:
   typename FiniteElement<1,spacedim>::InternalDataBase *
   get_data (const UpdateFlags,
             const Mapping<1,spacedim> &mapping,
-            const Quadrature<1> &quadrature) const;
+            const Quadrature<1> &quadrature) const
+  {
+    return new typename FiniteElement<1, spacedim>::InternalDataBase;
+  }
 
   typename FiniteElement<1,spacedim>::InternalDataBase *
-  get_face_data (const UpdateFlags,
-                 const Mapping<1,spacedim> &mapping,
-                 const Quadrature<0> &quadrature) const;
+  get_face_data(const UpdateFlags update_flags,
+                const Mapping<1,spacedim> &mapping,
+                const Quadrature<0> &quadrature) const
+  {
+    // generate a new data object and initialize some fields
+    typename FiniteElement<1, spacedim>::InternalDataBase *data =
+      new typename FiniteElement<1, spacedim>::InternalDataBase;
+
+    // check what needs to be initialized only once and what on every
+    // cell/face/subface we visit
+    data->update_once = update_once(update_flags);
+    data->update_each = update_each(update_flags);
+    data->update_flags = data->update_once | data->update_each;
+
+    const UpdateFlags flags(data->update_flags);
+    const unsigned int n_q_points = quadrature.size();
+    AssertDimension(n_q_points, 1);
+    (void)n_q_points;
+
+    // No derivatives of this element are implemented.
+    if (flags & update_gradients || flags & update_hessians)
+      {
+        Assert(false, ExcNotImplemented());
+      }
+
+    return data;
+  }
 
   typename FiniteElement<1,spacedim>::InternalDataBase *
-  get_subface_data (const UpdateFlags,
-                    const Mapping<1,spacedim> &mapping,
-                    const Quadrature<0> &quadrature) const;
+  get_subface_data(const UpdateFlags update_flags,
+                   const Mapping<1,spacedim> &mapping,
+                   const Quadrature<0> &quadrature) const
+  {
+    return get_face_data(update_flags, mapping, quadrature);
+  }
 
   virtual
   void
index 41ba86416adc3c4258511ca7bbffa7ff4de4461a..fb12ed105a19decf6a57421ce6b140657bdf35d4 100644 (file)
@@ -18,6 +18,7 @@
 
 
 #include <deal.II/fe/fe.h>
+#include <deal.II/base/quadrature.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -63,6 +64,7 @@ DEAL_II_NAMESPACE_OPEN
  *
  * @author Ralf Hartmann 2004, Guido Kanschat, 2009
  */
+
 template <class POLY, int dim=POLY::dimension, int spacedim=dim>
 class FE_Poly : public FiniteElement<dim,spacedim>
 {
@@ -164,9 +166,77 @@ protected:
 
   virtual
   typename FiniteElement<dim,spacedim>::InternalDataBase *
-  get_data (const UpdateFlags,
-            const Mapping<dim,spacedim> &mapping,
-            const Quadrature<dim> &quadrature) const;
+  get_data(const UpdateFlags update_flags,
+           const Mapping<dim,spacedim> &mapping,
+           const Quadrature<dim> &quadrature) const
+  {
+    // generate a new data object and
+    // initialize some fields
+    InternalData *data = new InternalData;
+
+    // check what needs to be
+    // initialized only once and what
+    // on every cell/face/subface we
+    // visit
+    data->update_once = update_once(update_flags);
+    data->update_each = update_each(update_flags);
+    data->update_flags = data->update_once | data->update_each;
+
+    const UpdateFlags flags(data->update_flags);
+    const unsigned int n_q_points = quadrature.size();
+
+    // some scratch arrays
+    std::vector<double> values(0);
+    std::vector<Tensor<1, dim> > grads(0);
+    std::vector<Tensor<2, dim> > grad_grads(0);
+
+    // initialize fields only if really
+    // necessary. otherwise, don't
+    // allocate memory
+    if (flags & update_values)
+      {
+        values.resize(this->dofs_per_cell);
+        data->shape_values.resize(this->dofs_per_cell,
+                                  std::vector<double>(n_q_points));
+      }
+
+    if (flags & update_gradients)
+      {
+        grads.resize(this->dofs_per_cell);
+        data->shape_gradients.resize(this->dofs_per_cell,
+                                     std::vector<Tensor<1, dim> >(n_q_points));
+      }
+
+    // if second derivatives through
+    // finite differencing is required,
+    // then initialize some objects for
+    // that
+    if (flags & update_hessians)
+      data->initialize_2nd(this, mapping, quadrature);
+
+    // next already fill those fields
+    // of which we have information by
+    // now. note that the shape
+    // gradients are only those on the
+    // unit cell, and need to be
+    // transformed when visiting an
+    // actual cell
+    if (flags & (update_values | update_gradients))
+      for (unsigned int i = 0; i<n_q_points; ++i)
+        {
+          poly_space.compute(quadrature.point(i),
+                             values, grads, grad_grads);
+
+          if (flags & update_values)
+            for (unsigned int k = 0; k<this->dofs_per_cell; ++k)
+              data->shape_values[k][i] = values[k];
+
+          if (flags & update_gradients)
+            for (unsigned int k = 0; k<this->dofs_per_cell; ++k)
+              data->shape_gradients[k][i] = grads[k];
+        }
+    return data;
+  }
 
   virtual
   void
index ce2a7bf4d66849f6cace81ab9f673d804d7f7b81..63d018447427b09d9de870206b21d3b80899a495 100644 (file)
@@ -161,96 +161,6 @@ FE_Poly<POLY,dim,spacedim>::update_each (const UpdateFlags flags) const
 
 
 
-//---------------------------------------------------------------------------
-// Data field initialization
-//---------------------------------------------------------------------------
-
-template <class POLY, int dim, int spacedim>
-typename FiniteElement<dim,spacedim>::InternalDataBase *
-FE_Poly<POLY,dim,spacedim>::get_data (const UpdateFlags      update_flags,
-                                      const Mapping<dim,spacedim> &,
-                                      const Quadrature<dim> &quadrature) const
-{
-  // generate a new data object and
-  // initialize some fields
-  InternalData *data = new InternalData;
-
-  // check what needs to be
-  // initialized only once and what
-  // on every cell/face/subface we
-  // visit
-  data->update_once = update_once(update_flags);
-  data->update_each = update_each(update_flags);
-  data->update_flags = data->update_once | data->update_each;
-
-  const UpdateFlags flags(data->update_flags);
-  const unsigned int n_q_points = quadrature.size();
-
-  // some scratch arrays
-  std::vector<double> values(0);
-  std::vector<Tensor<1,dim> > grads(0);
-  std::vector<Tensor<2,dim> > grad_grads(0);
-
-  // initialize fields only if really
-  // necessary. otherwise, don't
-  // allocate memory
-  if (flags & update_values)
-    {
-      values.resize (this->dofs_per_cell);
-      data->shape_values.resize (this->dofs_per_cell,
-                                 std::vector<double> (n_q_points));
-    }
-
-  if (flags & update_gradients)
-    {
-      grads.resize (this->dofs_per_cell);
-      data->shape_gradients.resize (this->dofs_per_cell,
-                                    std::vector<Tensor<1,dim> > (n_q_points));
-    }
-
-  // if second derivatives through
-  // finite differencing is required,
-  // then initialize some objects for
-  // that
-  if (flags & update_hessians)
-    {
-      grad_grads.resize (this->dofs_per_cell);
-      data->shape_hessians.resize (this->dofs_per_cell,
-                                   std::vector<Tensor<2,dim> > (n_q_points));
-      data->untransformed_shape_hessians.resize (n_q_points);
-    }
-
-  // next already fill those fields
-  // of which we have information by
-  // now. note that the shape
-  // gradients are only those on the
-  // unit cell, and need to be
-  // transformed when visiting an
-  // actual cell
-  if (flags & (update_values | update_gradients | update_hessians))
-    for (unsigned int i=0; i<n_q_points; ++i)
-      {
-        poly_space.compute(quadrature.point(i),
-                           values, grads, grad_grads);
-
-        if (flags & update_values)
-          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
-            data->shape_values[k][i] = values[k];
-
-        if (flags & update_gradients)
-          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
-            data->shape_gradients[k][i] = grads[k];
-
-        if (flags & update_hessians)
-          for (unsigned int k=0; k<this->dofs_per_cell; ++k)
-            data->shape_hessians[k][i] = grad_grads[k];
-      }
-  return data;
-}
-
-
-
-
 //---------------------------------------------------------------------------
 // Fill data of FEValues
 //---------------------------------------------------------------------------
index e0997eb93213843da09a9d72edb04a504877ba48..628a2d4fff3c169e3643fb20a71693ceb14cddcc 100644 (file)
@@ -78,17 +78,72 @@ protected:
   typename FiniteElement<dim,spacedim>::InternalDataBase *
   get_data (const UpdateFlags,
             const Mapping<dim,spacedim> &mapping,
-            const Quadrature<dim> &quadrature) const;
+            const Quadrature<dim> &quadrature) const
+  {
+    InternalData *data = new InternalData;
+    return data;
+  }
 
   typename FiniteElement<dim,spacedim>::InternalDataBase *
-  get_face_data (const UpdateFlags,
-                 const Mapping<dim,spacedim> &mapping,
-                 const Quadrature<dim-1>& quadrature) const;
+  get_face_data(const UpdateFlags update_flags,
+                const Mapping<dim,spacedim> &mapping,
+                const Quadrature<dim-1>& quadrature) const
+  {
+    // generate a new data object and
+    // initialize some fields
+    InternalData *data = new InternalData;
+
+    // check what needs to be
+    // initialized only once and what
+    // on every cell/face/subface we
+    // visit
+    data->update_once = update_once(update_flags);
+    data->update_each = update_each(update_flags);
+    data->update_flags = data->update_once | data->update_each;
+
+    const UpdateFlags flags(data->update_flags);
+    const unsigned int n_q_points = quadrature.size();
+
+    // some scratch arrays
+    std::vector<double> values(0);
+    std::vector<Tensor<1, dim - 1> > grads(0);
+    std::vector<Tensor<2, dim - 1> > grad_grads(0);
+
+    // initialize fields only if really
+    // necessary. otherwise, don't
+    // allocate memory
+    if (flags & update_values)
+      {
+        values.resize(poly_space.n());
+        data->shape_values.resize(poly_space.n(),
+                                  std::vector<double>(n_q_points));
+        for (unsigned int i = 0; i<n_q_points; ++i)
+          {
+            poly_space.compute(quadrature.point(i),
+                               values, grads, grad_grads);
+
+            for (unsigned int k = 0; k<poly_space.n(); ++k)
+              data->shape_values[k][i] = values[k];
+          }
+      }
+    // No derivatives of this element
+    // are implemented.
+    if (flags & update_gradients || flags & update_hessians)
+      {
+        Assert(false, ExcNotImplemented());
+      }
+
+    return data;
+  }
 
   typename FiniteElement<dim,spacedim>::InternalDataBase *
-  get_subface_data (const UpdateFlags,
-                    const Mapping<dim,spacedim> &mapping,
-                    const Quadrature<dim-1>& quadrature) const;
+  get_subface_data(const UpdateFlags update_flags,
+                   const Mapping<dim,spacedim> &mapping,
+                   const Quadrature<dim-1>& quadrature) const
+  {
+    return get_face_data(update_flags, mapping,
+                         QProjector<dim - 1>::project_to_all_children(quadrature));
+  }
 
   virtual
   void
index ca0a4037b7225148cd8de55163f505de3ec93a71..786f0bd3e935d3c683df135b53d79d4cdb424d65 100644 (file)
@@ -78,94 +78,6 @@ FE_PolyFace<POLY,dim,spacedim>::update_each (const UpdateFlags flags) const
 }
 
 
-
-//---------------------------------------------------------------------------
-// Data field initialization
-//---------------------------------------------------------------------------
-
-template <class POLY, int dim, int spacedim>
-typename FiniteElement<dim,spacedim>::InternalDataBase *
-FE_PolyFace<POLY,dim,spacedim>::get_data (
-  const UpdateFlags,
-  const Mapping<dim,spacedim> &,
-  const Quadrature<dim> &) const
-{
-  InternalData *data = new InternalData;
-  return data;
-}
-
-
-template <class POLY, int dim, int spacedim>
-typename FiniteElement<dim,spacedim>::InternalDataBase *
-FE_PolyFace<POLY,dim,spacedim>::get_face_data (
-  const UpdateFlags update_flags,
-  const Mapping<dim,spacedim> &,
-  const Quadrature<dim-1>& quadrature) const
-{
-  // generate a new data object and
-  // initialize some fields
-  InternalData *data = new InternalData;
-
-  // check what needs to be
-  // initialized only once and what
-  // on every cell/face/subface we
-  // visit
-  data->update_once = update_once(update_flags);
-  data->update_each = update_each(update_flags);
-  data->update_flags = data->update_once | data->update_each;
-
-  const UpdateFlags flags(data->update_flags);
-  const unsigned int n_q_points = quadrature.size();
-
-  // some scratch arrays
-  std::vector<double> values(0);
-  std::vector<Tensor<1,dim-1> > grads(0);
-  std::vector<Tensor<2,dim-1> > grad_grads(0);
-
-  // initialize fields only if really
-  // necessary. otherwise, don't
-  // allocate memory
-  if (flags & update_values)
-    {
-      values.resize (poly_space.n());
-      data->shape_values.resize (poly_space.n(),
-                                 std::vector<double> (n_q_points));
-      for (unsigned int i=0; i<n_q_points; ++i)
-        {
-          poly_space.compute(quadrature.point(i),
-                             values, grads, grad_grads);
-
-          for (unsigned int k=0; k<poly_space.n(); ++k)
-            data->shape_values[k][i] = values[k];
-        }
-    }
-  // No derivatives of this element
-  // are implemented.
-  if (flags & update_gradients || flags & update_hessians)
-    {
-      Assert(false, ExcNotImplemented());
-    }
-
-  return data;
-}
-
-
-
-template <class POLY, int dim, int spacedim>
-typename FiniteElement<dim,spacedim>::InternalDataBase *
-FE_PolyFace<POLY,dim,spacedim>::get_subface_data (
-  const UpdateFlags flags,
-  const Mapping<dim,spacedim> &mapping,
-  const Quadrature<dim-1>& quadrature) const
-{
-  return get_face_data (flags, mapping,
-                        QProjector<dim-1>::project_to_all_children(quadrature));
-}
-
-
-
-
-
 //---------------------------------------------------------------------------
 // Fill data of FEValues
 //---------------------------------------------------------------------------
index 0c9c5ed9b394ebdd296fd36ac1853a9b571e21d9..63b2e4b6e5f7ef9468f2d43a55d57f226f78c404 100644 (file)
@@ -177,9 +177,104 @@ protected:
 
   virtual
   typename FiniteElement<dim,spacedim>::InternalDataBase *
-  get_data (const UpdateFlags,
-            const Mapping<dim,spacedim> &mapping,
-            const Quadrature<dim> &quadrature) const;
+  get_data(const UpdateFlags update_flags,
+           const Mapping<dim,spacedim> &mapping,
+           const Quadrature<dim> &quadrature) const
+  {
+    // generate a new data object and
+    // initialize some fields
+    InternalData *data = new InternalData;
+
+    // check what needs to be
+    // initialized only once and what
+    // on every cell/face/subface we
+    // visit
+    data->update_once = update_once(update_flags);
+    data->update_each = update_each(update_flags);
+    data->update_flags = data->update_once | data->update_each;
+
+    const UpdateFlags flags(data->update_flags);
+    const unsigned int n_q_points = quadrature.size();
+
+    // some scratch arrays
+    std::vector<Tensor<1, dim> > values(0);
+    std::vector<Tensor<2, dim> > grads(0);
+    std::vector<Tensor<3, dim> > grad_grads(0);
+
+    // initialize fields only if really
+    // necessary. otherwise, don't
+    // allocate memory
+    if (flags & update_values)
+      {
+        values.resize(this->dofs_per_cell);
+        data->shape_values.resize(this->dofs_per_cell);
+        for (unsigned int i = 0; i<this->dofs_per_cell; ++i)
+          data->shape_values[i].resize(n_q_points);
+      }
+
+    if (flags & update_gradients)
+      {
+        grads.resize(this->dofs_per_cell);
+        data->shape_grads.resize(this->dofs_per_cell);
+        for (unsigned int i = 0; i<this->dofs_per_cell; ++i)
+          data->shape_grads[i].resize(n_q_points);
+      }
+
+    // if second derivatives through
+    // finite differencing is required,
+    // then initialize some objects for
+    // that
+    if (flags & update_hessians)
+      {
+        //      grad_grads.resize (this->dofs_per_cell);
+        data->initialize_2nd(this, mapping, quadrature);
+      }
+
+    // Compute shape function values
+    // and derivatives on the reference
+    // cell. Make sure, that for the
+    // node values N_i holds
+    // N_i(v_j)=\delta_ij for all basis
+    // functions v_j
+    if (flags & (update_values | update_gradients))
+      for (unsigned int k = 0; k<n_q_points; ++k)
+        {
+          poly_space.compute(quadrature.point(k),
+                             values, grads, grad_grads);
+
+          if (flags & update_values)
+            {
+              if (inverse_node_matrix.n_cols() == 0)
+                for (unsigned int i = 0; i<this->dofs_per_cell; ++i)
+                  data->shape_values[i][k] = values[i];
+              else
+                for (unsigned int i = 0; i<this->dofs_per_cell; ++i)
+                  {
+                    Tensor<1, dim> add_values;
+                    for (unsigned int j = 0; j<this->dofs_per_cell; ++j)
+                      add_values += inverse_node_matrix(j, i) * values[j];
+                    data->shape_values[i][k] = add_values;
+                  }
+            }
+
+          if (flags & update_gradients)
+            {
+              if (inverse_node_matrix.n_cols() == 0)
+                for (unsigned int i = 0; i<this->dofs_per_cell; ++i)
+                  data->shape_grads[i][k] = grads[i];
+              else
+                for (unsigned int i = 0; i<this->dofs_per_cell; ++i)
+                  {
+                    Tensor<2, dim> add_grads;
+                    for (unsigned int j = 0; j<this->dofs_per_cell; ++j)
+                      add_grads += inverse_node_matrix(j, i) * grads[j];
+                    data->shape_grads[i][k] = add_grads;
+                  }
+            }
+
+        }
+    return data;
+  }
 
   virtual
   void
index 68418bfc68b61cd4fb11835f3554487c09cb29f9..b3d9e2db48479417a6c10ffdf38f54d915b3438f 100644 (file)
@@ -458,63 +458,6 @@ FE_FaceQ<1,spacedim>::update_each (const UpdateFlags flags) const
 }
 
 
-
-template <int spacedim>
-typename FiniteElement<1,spacedim>::InternalDataBase *
-FE_FaceQ<1,spacedim>::get_data (
-  const UpdateFlags,
-  const Mapping<1,spacedim> &,
-  const Quadrature<1> &) const
-{
-  return new typename FiniteElement<1,spacedim>::InternalDataBase;
-}
-
-
-template <int spacedim>
-typename FiniteElement<1,spacedim>::InternalDataBase *
-FE_FaceQ<1,spacedim>::get_face_data (
-  const UpdateFlags update_flags,
-  const Mapping<1,spacedim> &,
-  const Quadrature<0> &quadrature) const
-{
-  // generate a new data object and initialize some fields
-  typename FiniteElement<1,spacedim>::InternalDataBase *data =
-    new typename FiniteElement<1,spacedim>::InternalDataBase;
-
-  // check what needs to be initialized only once and what on every
-  // cell/face/subface we visit
-  data->update_once = update_once(update_flags);
-  data->update_each = update_each(update_flags);
-  data->update_flags = data->update_once | data->update_each;
-
-  const UpdateFlags flags(data->update_flags);
-  const unsigned int n_q_points = quadrature.size();
-  AssertDimension(n_q_points, 1);
-  (void)n_q_points;
-
-  // No derivatives of this element are implemented.
-  if (flags & update_gradients || flags & update_hessians)
-    {
-      Assert(false, ExcNotImplemented());
-    }
-
-  return data;
-}
-
-
-
-template <int spacedim>
-typename FiniteElement<1,spacedim>::InternalDataBase *
-FE_FaceQ<1,spacedim>::get_subface_data (
-  const UpdateFlags flags,
-  const Mapping<1,spacedim> &mapping,
-  const Quadrature<0> &quadrature) const
-{
-  return get_face_data (flags, mapping, quadrature);
-}
-
-
-
 template <int spacedim>
 void
 FE_FaceQ<1,spacedim>::
index 57bdd8b80be1df977db7d02d6b5b4d99f30da9ab..fc9e1e752ac5be3f42c40e8ba45b2c159db2bcc3 100644 (file)
@@ -273,127 +273,6 @@ FE_PolyTensor<POLY,dim,spacedim>::shape_grad_grad_component (const unsigned int
 }
 
 
-
-//---------------------------------------------------------------------------
-// Data field initialization
-//---------------------------------------------------------------------------
-
-template <class POLY, int dim, int spacedim>
-typename FiniteElement<dim,spacedim>::InternalDataBase *
-FE_PolyTensor<POLY,dim,spacedim>::get_data (
-  const UpdateFlags      update_flags,
-  const Mapping<dim,spacedim>    &mapping,
-  const Quadrature<dim> &quadrature) const
-{
-  // generate a new data object and
-  // initialize some fields
-  InternalData *data = new InternalData;
-
-  // check what needs to be
-  // initialized only once and what
-  // on every cell/face/subface we
-  // visit
-  data->update_once = update_once(update_flags);
-  data->update_each = update_each(update_flags);
-  data->update_flags = data->update_once | data->update_each;
-
-  const UpdateFlags flags(data->update_flags);
-  const unsigned int n_q_points = quadrature.size();
-
-  // some scratch arrays
-  std::vector<Tensor<1,dim> > values(0);
-  std::vector<Tensor<2,dim> > grads(0);
-  std::vector<Tensor<3,dim> > grad_grads(0);
-
-  data->sign_change.resize (this->dofs_per_cell);
-
-  // initialize fields only if really
-  // necessary. otherwise, don't
-  // allocate memory
-  if (flags & update_values)
-    {
-      values.resize (this->dofs_per_cell);
-      data->shape_values.resize (this->dofs_per_cell);
-      for (unsigned int i=0; i<this->dofs_per_cell; ++i)
-        data->shape_values[i].resize (n_q_points);
-      if (mapping_type != mapping_none)
-        data->transformed_shape_values.resize (n_q_points);
-    }
-
-  if (flags & update_gradients)
-    {
-      grads.resize (this->dofs_per_cell);
-      data->shape_grads.resize (this->dofs_per_cell);
-      for (unsigned int i=0; i<this->dofs_per_cell; ++i)
-        data->shape_grads[i].resize (n_q_points);
-      data->transformed_shape_grads.resize (n_q_points);
-      if ( (mapping_type == mapping_raviart_thomas)
-           ||
-           (mapping_type == mapping_piola)
-           ||
-           (mapping_type == mapping_nedelec))
-        data->untransformed_shape_grads.resize(n_q_points);
-    }
-
-  // if second derivatives through
-  // finite differencing is required,
-  // then initialize some objects for
-  // that
-  if (flags & update_hessians)
-    {
-//      grad_grads.resize (this->dofs_per_cell);
-      data->initialize_2nd (this, mapping, quadrature);
-    }
-
-  // Compute shape function values
-  // and derivatives on the reference
-  // cell. Make sure, that for the
-  // node values N_i holds
-  // N_i(v_j)=\delta_ij for all basis
-  // functions v_j
-  if (flags & (update_values | update_gradients))
-    for (unsigned int k=0; k<n_q_points; ++k)
-      {
-        poly_space.compute(quadrature.point(k),
-                           values, grads, grad_grads);
-
-        if (flags & update_values)
-          {
-            if (inverse_node_matrix.n_cols() == 0)
-              for (unsigned int i=0; i<this->dofs_per_cell; ++i)
-                data->shape_values[i][k] = values[i];
-            else
-              for (unsigned int i=0; i<this->dofs_per_cell; ++i)
-                {
-                  Tensor<1,dim> add_values;
-                  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
-                    add_values += inverse_node_matrix(j,i) * values[j];
-                  data->shape_values[i][k] = add_values;
-                }
-          }
-
-        if (flags & update_gradients)
-          {
-            if (inverse_node_matrix.n_cols() == 0)
-              for (unsigned int i=0; i<this->dofs_per_cell; ++i)
-                data->shape_grads[i][k] = grads[i];
-            else
-              for (unsigned int i=0; i<this->dofs_per_cell; ++i)
-                {
-                  Tensor<2,dim> add_grads;
-                  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
-                    add_grads += inverse_node_matrix(j,i) * grads[j];
-                  data->shape_grads[i][k] = add_grads;
-                }
-          }
-
-      }
-  return data;
-}
-
-
-
-
 //---------------------------------------------------------------------------
 // Fill data of FEValues
 //---------------------------------------------------------------------------

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.