]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Removed a couple of functions.
authorLuca Heltai <luca.heltai@sissa.it>
Thu, 3 Mar 2016 09:46:32 +0000 (10:46 +0100)
committerLuca Heltai <luca.heltai@sissa.it>
Wed, 6 Apr 2016 11:00:14 +0000 (13:00 +0200)
include/deal.II/fe/mapping_manifold.h
source/fe/mapping_manifold.cc
tests/fe/mapping_manifold_02.cc

index 9d2c5b111a4633180f7f14c6b4a039c6346f470c..b7b7b26dfedc296da5ded0f225fd23426c853cc6 100644 (file)
@@ -205,18 +205,17 @@ public:
     void
     compute_manifold_quadrature_weights(const Quadrature<dim> &quadrature);
 
-
     /**
      * Return an estimate (in bytes) or the memory consumption of this object.
      */
     virtual std::size_t memory_consumption () const;
 
     /**
-     * Values of manifold quadrature formulas.
+     * Store the current cell.
      *
      * Computed each.
      */
-    std::vector<Quadrature<spacedim> > cell_manifold_quadratures;
+    typename Triangulation<dim,spacedim>::cell_iterator current_cell;
 
     /**
      * Values of quadrature weights for manifold quadrature formulas.
index 22ec5554bca215658bfda0a6d24c9e7e847bb6b0..2ae478fe977a76c23c83ba68aae28156aac7a4d5 100644 (file)
@@ -78,13 +78,10 @@ initialize (const UpdateFlags      update_flags,
 
   const unsigned int n_q_points = q.size();
 
-  // Update the weights used in the Manifold Quadrature formulas
-  compute_manifold_quadrature_weights(q);
-
   // see if we need the (transformation) shape function values
   // and/or gradients and resize the necessary arrays
   if (this->update_each & update_quadrature_points)
-    cell_manifold_quadratures.resize(q.size());
+    compute_manifold_quadrature_weights(q);
 
   // if (this->update_each & (update_covariant_transformation
   //                          | update_contravariant_transformation
@@ -127,7 +124,6 @@ initialize (const UpdateFlags      update_flags,
 }
 
 
-
 template <int dim, int spacedim>
 void
 MappingManifold<dim,spacedim>::InternalData::
@@ -377,23 +373,26 @@ namespace internal
      */
     template <int dim, int spacedim>
     void
-    maybe_compute_q_points (const typename QProjector<dim>::DataSetDescriptor                 data_set,
+    maybe_compute_q_points (const typename dealii::Triangulation<dim,spacedim>::cell_iterator &cell,
+                            const typename QProjector<dim>::DataSetDescriptor                 data_set,
                             const typename dealii::MappingManifold<dim,spacedim>::InternalData      &data,
                             std::vector<Point<spacedim> >                                     &quadrature_points)
     {
       const UpdateFlags update_flags = data.update_each;
 
+      AssertDimension(quadrature_points.size(),
+                      data.cell_manifold_quadratures_weights.size());
+
+      std::vector<Point<spacedim> > vertices;
+      for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
+        vertices[v] = cell->vertex(v);
+
       if (update_flags & update_quadrature_points)
         {
           for (unsigned int point=0; point<quadrature_points.size(); ++point)
             {
-              const double *shape = &data.shape(point+data_set,0);
-              Point<spacedim> result = (shape[0] *
-                                        data.mapping_support_points[0]);
-              for (unsigned int k=1; k<data.n_shape_functions; ++k)
-                for (unsigned int i=0; i<spacedim; ++i)
-                  result[i] += shape[k] * data.mapping_support_points[k][i];
-              quadrature_points[point] = result;
+              quadrature_points[point] = cell->get_manifold().
+                                         get_new_point(Quadrature<spacedim>(vertices, data.cell_manifold_quadratures_weights[point]));
             }
         }
     }
@@ -919,31 +918,18 @@ fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                 const typename Mapping<dim,spacedim>::InternalDataBase    &internal_data,
                 internal::FEValues::MappingRelatedData<dim,spacedim>      &output_data) const
 {
-//   // ensure that the following static_cast is really correct:
-//   Assert (dynamic_cast<const InternalData *>(&internal_data) != 0,
-//           ExcInternalError());
-//   const InternalData &data = static_cast<const InternalData &>(internal_data);
+  // ensure that the following static_cast is really correct:
+  Assert (dynamic_cast<const InternalData *>(&internal_data) != 0,
+          ExcInternalError());
+  const InternalData &data = static_cast<const InternalData &>(internal_data);
 
-//   const unsigned int n_q_points=quadrature.size();
+  const unsigned int n_q_points=quadrature.size();
 
-//   // if necessary, recompute the support points of the transformation of this cell
-//   // (note that we need to first check the triangulation pointer, since otherwise
-//   // the second test might trigger an exception if the triangulations are not the
-//   // same)
-//   if ((data.mapping_support_points.size() == 0)
-//       ||
-//       (&cell->get_triangulation() !=
-//        &data.cell_of_current_support_points->get_triangulation())
-//       ||
-//       (cell != data.cell_of_current_support_points))
-//     {
-//       data.mapping_support_points = this->compute_mapping_support_points(cell);
-//       data.cell_of_current_support_points = cell;
-//     }
+  internal::maybe_compute_q_points<dim,spacedim> (cell,
+                                                  QProjector<dim>::DataSetDescriptor::cell (),
+                                                  data,
+                                                  output_data.quadrature_points);
 
-//   internal::maybe_compute_q_points<dim,spacedim> (QProjector<dim>::DataSetDescriptor::cell (),
-//                                                   data,
-//                                                   output_data.quadrature_points);
 //   internal::maybe_update_Jacobians<dim,spacedim> (cell_similarity,
 //                                                   QProjector<dim>::DataSetDescriptor::cell (),
 //                                                   data);
@@ -1258,7 +1244,8 @@ namespace internal
                             const typename dealii::MappingManifold<dim,spacedim>::InternalData      &data,
                             internal::FEValues::MappingRelatedData<dim,spacedim>              &output_data)
     {
-      maybe_compute_q_points<dim,spacedim> (data_set,
+      maybe_compute_q_points<dim,spacedim> (cell,
+                                            data_set,
                                             data,
                                             output_data.quadrature_points);
       maybe_update_Jacobians<dim,spacedim> (CellSimilarity::none,
index 2bb80f5bbaba6fadbec1d34cb0e61bdeaa36b2b0..4fc61bb5179e549be0c079bb6ec996abfe904298 100644 (file)
@@ -84,10 +84,11 @@ void test()
         {
           const Point<spacedim> pq = map_manifold.transform_unit_to_real_cell(cell,q_points_from_quadrature[q]);
 
-          if(pq.distance(q_points_from_fe_values[q]) > 1e-10) {
-           deallog << "Expected: " << pq << ", got: "
-                   << q_points_from_fe_values[q] << std::endl;
-         }
+          if (pq.distance(q_points_from_fe_values[q]) > 1e-10)
+            {
+              deallog << "Expected: " << pq << ", got: "
+                      << q_points_from_fe_values[q] << std::endl;
+            }
         }
     }
   deallog << "OK" << std::endl;

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.