]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Extend the Patch class to allow for different space and object dimensions, such as...
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 7 Sep 2000 13:00:33 +0000 (13:00 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 7 Sep 2000 13:00:33 +0000 (13:00 +0000)
git-svn-id: https://svn.dealii.org/trunk@3303 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/base/include/base/data_out_base.h
deal.II/base/source/data_out_base.all_dimensions.cc
deal.II/base/source/data_out_base.cc
deal.II/doc/news/2000/c-3-0.html

index db072abadb3e9d99c7efaea23e8bc5c3d75d8e70..da9b67032fce7e2e900c785a95002812cbbaca7b 100644 (file)
@@ -80,8 +80,24 @@ class ParameterHandler;
  * i.e. the z-coordinate runs fastest, then the y-coordinate, then x (if there
  * are that many space directions).
  *
- * The @p{Patch} class takes a template parameter denoting the space dimension
- * in which this patch operates.
+ *
+ * @sect3{Generalized patches}
+ *
+ * In general, the patches as explained above might be too
+ * restricted. For example, one might want to draw only the outer
+ * faces of a domain in a three-dimensional computation, if one is not
+ * interested in what happens inside. Then, the objects that should be
+ * drawn are two-dimensional in a three-dimensional world. The
+ * @p{Patch} class and associated output functions handle these
+ * cases. The @p{Patch} class therefore takes two template parameters,
+ * the first, named @p{dim} denoting the dimension of the object (in
+ * the above example, this would be two), while the second, named
+ * @p{spacedim}, denotes the dimension of the embedding space (this
+ * would be three). The corner points of a patch have the dimension of
+ * the space, while their number is determined by the dimension of the
+ * patch. By default, the second template parameter has the same value
+ * as the first, which would correspond to outputting a cell, rather
+ * than a face or something else.
  *
  *
  * @sect2{Supported output formats}
@@ -155,9 +171,8 @@ class ParameterHandler;
  * with values in the third direction taken from a data vector.
  *
  * The output uses two different povray-objects:
-
- * @begin{itemize}
  *
+ * @begin{itemize}
  * @item @p{BICUBIC_PATCH}
  * A @p{bicubic_patch} is a 3-dimensional Bezier patch. It consists of 16 Points
  * describing the surface. The 4 corner points are touched by the object,
@@ -192,6 +207,7 @@ class ParameterHandler;
  * If the external file "data.inc" is used, the path to this file has to be
  * included in the povray options.
  *
+ *
  * @sect3{EPS format}
  *
  * Output in this format circumvents the use of auxiliary graphic programs
@@ -267,7 +283,7 @@ class ParameterHandler;
  * details of the viewpoint, light source, etc.
  *
  *
- * @author Wolfgang Bangerth 1999; EPS output based on an earlier implementation by Stefan Nauber for the old DataOut class; Povray output by Thomas Richter 1999
+ * @author Wolfgang Bangerth 1999, 2000; EPS output based on an earlier implementation by Stefan Nauber for the old DataOut class; Povray output by Thomas Richter 1999
  */
 class DataOutBase 
 {
@@ -298,7 +314,7 @@ class DataOutBase
                                      * @end{verbatim}
                                      * @author Wolfgang Bangerth
                                      */
-    template <int dim>
+    template <int dim, int spacedim=dim>
     struct Patch
     {
                                         /**
@@ -311,7 +327,7 @@ class DataOutBase
                                          * is the same as for cells
                                          * in the triangulation.
                                          */
-       Point<dim> vertices[GeometryInfo<dim>::vertices_per_cell];
+       Point<spacedim> vertices[GeometryInfo<dim>::vertices_per_cell];
        
                                         /**
                                          * Number of subdivisions with
@@ -401,6 +417,15 @@ class DataOutBase
                                          * @p{n_subdivisions} to one.
                                          */
        Patch ();
+
+                                        /**
+                                         * Exception
+                                         */
+       DeclException2 (ExcInvalidCombinationOfDimensions,
+                       int, int,
+                       << "It is not possible to have a structural dimension of " << arg1
+                       << " to be larger than the space dimension of the surrounding"
+                       << " space " << arg2);
     };
 
 
@@ -869,11 +894,11 @@ class DataOutBase
                                      * documentation for more information
                                      * on the parameters.
                                      */
-    template <int dim>
-    static void write_ucd (const vector<Patch<dim> > &patches,
-                          const vector<string>      &data_names,
-                          const UcdFlags            &flags,
-                          ostream                   &out);
+    template <int dim, int spacedim>
+    static void write_ucd (const vector<Patch<dim,spacedim> > &patches,
+                          const vector<string>               &data_names,
+                          const UcdFlags                     &flags,
+                          ostream                            &out);
 
                                     /**
                                      * Write the given list of patches
@@ -882,11 +907,11 @@ class DataOutBase
                                      * documentation for more information
                                      * on the parameters.
                                      */
-    template <int dim>
-    static void write_gnuplot (const vector<Patch<dim> > &patches,
-                              const vector<string>      &data_names,
-                              const GnuplotFlags        &flags,
-                              ostream                   &out);
+    template <int dim, int spacedim>
+    static void write_gnuplot (const vector<Patch<dim,spacedim> > &patches,
+                              const vector<string>               &data_names,
+                              const GnuplotFlags                 &flags,
+                              ostream                            &out);
 
                                     /**
                                      * Write the given list of patches
@@ -895,11 +920,11 @@ class DataOutBase
                                      * documentation for more information
                                      * on the parameters.
                                      */
-    template <int dim>
-    static void write_povray (const vector<Patch<dim> > &patches,
-                             const vector<string>      &data_names,
-                             const PovrayFlags         &flags,
-                             ostream                   &out);
+    template <int dim, int spacedim>
+    static void write_povray (const vector<Patch<dim,spacedim> > &patches,
+                             const vector<string>               &data_names,
+                             const PovrayFlags                  &flags,
+                             ostream                            &out);
 
                                     /**
                                      * Write the given list of patches
@@ -908,11 +933,11 @@ class DataOutBase
                                      * documentation for more information
                                      * on the parameters.
                                      */
-    template <int dim>
-    static void write_eps (const vector<Patch<dim> > &patches,
-                          const vector<string>      &data_names,
-                          const EpsFlags            &flags,
-                          ostream                   &out);
+    template <int dim, int spacedim>
+    static void write_eps (const vector<Patch<dim,spacedim> > &patches,
+                          const vector<string>               &data_names,
+                          const EpsFlags                     &flags,
+                          ostream                            &out);
 
                                     /**
                                      * Write the given list of patches
@@ -921,11 +946,11 @@ class DataOutBase
                                      * documentation for more information
                                      * on the parameters.
                                      */
-    template <int dim>
-    static void write_gmv (const vector<Patch<dim> > &patches,
-                          const vector<string>      &data_names,
-                          const GmvFlags            &flags,
-                          ostream                   &out);
+    template <int dim, int spacedim>
+    static void write_gmv (const vector<Patch<dim,spacedim> > &patches,
+                          const vector<string>               &data_names,
+                          const GmvFlags                     &flags,
+                          ostream                            &out);
 
 
                                     /**
@@ -1023,12 +1048,15 @@ class DataOutBase
                                      * thus moved into this separate
                                      * function.
                                      */
-    template <int dim>
-    static void write_gmv_reorder_data_vectors (const vector<Patch<dim> > &patches,
-                                               vector<vector<double> >   &data_vectors);
+    template <int dim, int spacedim>
+    static void
+    write_gmv_reorder_data_vectors (const vector<Patch<dim,spacedim> > &patches,
+                                   vector<vector<double> >            &data_vectors);
 };
 
 
+
+
 /**
  * This class is the interface to the @p{DataOutBase} class, as already its name
  * might suggest. It does not offer much functionality apart from a way to
@@ -1141,7 +1169,7 @@ class DataOutBase
  *
  * @author Wolfgang Bangerth, 1999
  */
-template <int dim>
+template <int dim, int spacedim=dim>
 class DataOutInterface : private DataOutBase
 {
   public:
@@ -1342,7 +1370,8 @@ class DataOutInterface : private DataOutBase
                                      * allow the output functions to
                                      * know what they shall print.
                                      */
-    virtual const vector<DataOutBase::Patch<dim> > & get_patches () const = 0;
+    virtual const vector<DataOutBase::Patch<dim,spacedim> > &
+    get_patches () const = 0;
 
                                     /**
                                      * Abstract virtual function through
index 9c85ace5eeb0451eca05a8375a70630dff1babcd..2fcad5f315cf89d170434d57bb96c220bb46b796 100644 (file)
@@ -21,6 +21,7 @@ DataOutBase::UcdFlags::UcdFlags (const bool write_preamble) :
 {};
 
 
+
 DataOutBase::PovrayFlags::PovrayFlags (const bool smooth,
                                       const bool bicubic_patch,
                                       const bool external_data) :
@@ -30,26 +31,29 @@ DataOutBase::PovrayFlags::PovrayFlags (const bool smooth,
 {};
 
 
+
 void DataOutBase::UcdFlags::declare_parameters (ParameterHandler &prm)
 {
   prm.declare_entry ("Write preamble", "true", Patterns::Bool());
 };
 
 
+
 void DataOutBase::UcdFlags::parse_parameters (ParameterHandler &prm)
 {
   write_preamble = prm.get_bool ("Write preamble");
 };
 
 
+
 void DataOutBase::GnuplotFlags::declare_parameters (ParameterHandler &/*prm*/)
-{
-};
+{};
+
 
 
 void DataOutBase::GnuplotFlags::parse_parameters (ParameterHandler &/*prm*/)
-{
-};
+{};
+
 
 
 void DataOutBase::PovrayFlags::declare_parameters (ParameterHandler &prm)
@@ -63,6 +67,7 @@ void DataOutBase::PovrayFlags::declare_parameters (ParameterHandler &prm)
 };
 
 
+
 void DataOutBase::PovrayFlags::parse_parameters (ParameterHandler &prm)
 {
   smooth        = prm.get_bool ("Use smooth triangles");
@@ -71,6 +76,7 @@ void DataOutBase::PovrayFlags::parse_parameters (ParameterHandler &prm)
 };
 
 
+
 DataOutBase::EpsFlags::EpsFlags (const unsigned int  height_vector,
                                 const unsigned int  color_vector,
                                 const SizeType      size_type,
@@ -98,6 +104,7 @@ DataOutBase::EpsFlags::EpsFlags (const unsigned int  height_vector,
 {};
 
 
+
 DataOutBase::EpsFlags::RgbValues
 DataOutBase::EpsFlags::default_color_function (const double x,
                                               const double xmin,
@@ -196,6 +203,7 @@ DataOutBase::EpsFlags::grey_scale_color_function (const double x,
 };
 
 
+
 bool DataOutBase::EpsCell2d::operator < (const EpsCell2d &e) const
 {
                                   // note the "wrong" order in
@@ -204,6 +212,7 @@ bool DataOutBase::EpsCell2d::operator < (const EpsCell2d &e) const
 };
 
 
+
 void DataOutBase::EpsFlags::declare_parameters (ParameterHandler &prm)
 {
   prm.declare_entry ("Index of vector for height", "0",
@@ -233,6 +242,7 @@ void DataOutBase::EpsFlags::declare_parameters (ParameterHandler &prm)
 };
 
 
+
 void DataOutBase::EpsFlags::parse_parameters (ParameterHandler &prm)
 {
   height_vector = prm.get_integer ("Index of vector for height");
@@ -256,11 +266,11 @@ void DataOutBase::EpsFlags::parse_parameters (ParameterHandler &prm)
 };
 
 
+
 void DataOutBase::GmvFlags::declare_parameters (ParameterHandler &/*prm*/)
-{
-};
+{};
+
 
 
 void DataOutBase::GmvFlags::parse_parameters (ParameterHandler &/*prm*/)
-{
-};
+{};
index 879ca50f580b237041218c413d3b59b28a2a179f..4094f0ec6a9a0cba52e69430c5db7b0b07623703 100644 (file)
 #include <set>
 
 
-template <int dim>
-DataOutBase::Patch<dim>::Patch () :
+template <int dim, int spacedim>
+DataOutBase::Patch<dim,spacedim>::Patch () :
                n_subdivisions (1)
                                 // all the other data has a
                                 // constructor of its own
-{};
+{
+  Assert (dim<=spacedim, ExcInvalidCombinationOfDimensions(dim,spacedim));
+  Assert (spacedim<=3, ExcNotImplemented());
+};
 
 
-template <int dim>
-void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
-                            const vector<string>      &data_names,
-                            const UcdFlags            &flags,
-                            ostream                   &out) 
+
+template <int dim, int spacedim>
+void DataOutBase::write_ucd (const vector<Patch<dim,spacedim> > &patches,
+                            const vector<string>               &data_names,
+                            const UcdFlags                     &flags,
+                            ostream                            &out) 
 {
   AssertThrow (out, ExcIO());
 
@@ -46,7 +50,7 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
                                   // and cells for later use
   unsigned int n_cells = 0,
               n_nodes = 0;
-  for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+  for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
        patch!=patches.end(); ++patch)
     switch (dim)
       {
@@ -113,7 +117,7 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
     {
       unsigned int present_node = 1;
       
-      for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+      for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
           patch!=patches.end(); ++patch)
        {
          const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -125,11 +129,23 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
              case 1:
              {
                for (unsigned int i=0; i<n_subdivisions+1; ++i, ++present_node)
-                 out << present_node
-                     << "   "
-                     << ((patch->vertices[1](0) * i / n_subdivisions) +
-                         (patch->vertices[0](0) * (n_subdivisions-i) / n_subdivisions))
-                     << " 0 0\n";                        // fill with zeroes
+                 {
+                   out << present_node
+                       << "   ";
+
+                   const Point<spacedim>
+                     node = ((patch->vertices[1] * i / n_subdivisions) +
+                             (patch->vertices[0] * (n_subdivisions-i) / n_subdivisions));
+
+                                                    // write out coordinates
+                   for (unsigned int i=0; i<spacedim; ++i)
+                     out << node(i) << ' ';
+                                                    // fill with zeroes
+                   for (unsigned int i=spacedim; i<3; ++i)
+                     out << "0 ";
+                   out << endl;
+                 };
+               
                break;
              };
               
@@ -141,15 +157,24 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
                      const double x_frac = i * 1./n_subdivisions,
                                   y_frac = j * 1./n_subdivisions;
                      
+                     out << present_node
+                         << "   ";
+                     
                                                       // compute coordinates for
                                                       // this patch point
-                     out << present_node
-                         << "  "
-                         << (((patch->vertices[1] * x_frac) +
-                              (patch->vertices[0] * (1-x_frac))) * (1-y_frac) +
-                             ((patch->vertices[2] * x_frac) +
-                              (patch->vertices[3] * (1-x_frac))) * y_frac)
-                         << " 0\n";                   // fill with zeroes
+                     const Point<spacedim>
+                       node = (((patch->vertices[1] * x_frac) +
+                                (patch->vertices[0] * (1-x_frac))) * (1-y_frac) +
+                               ((patch->vertices[2] * x_frac) +
+                                (patch->vertices[3] * (1-x_frac))) * y_frac);
+                     
+                                                      // write out coordinates
+                     for (unsigned int i=0; i<spacedim; ++i)
+                       out << node(i) << ' ';
+                                                      // fill with zeroes
+                     for (unsigned int i=spacedim; i<3; ++i)
+                       out << "0 ";
+                     out << endl;
 
                      ++present_node;
                    };
@@ -177,16 +202,27 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
                                                         // compute coordinates for
                                                         // this patch point
                        out << present_node
-                           << "  "
-                           << ((((patch->vertices[1] * x_frac) +
-                                 (patch->vertices[0] * (1-x_frac))) * (1-y_frac) +
-                                ((patch->vertices[2] * x_frac) +
-                                 (patch->vertices[3] * (1-x_frac))) * y_frac)   * (1-z_frac) +
-                               (((patch->vertices[5] * x_frac) +
-                                 (patch->vertices[4] * (1-x_frac))) * (1-y_frac) +
-                                ((patch->vertices[6] * x_frac) +
-                                 (patch->vertices[7] * (1-x_frac))) * y_frac)   * z_frac)
-                           << endl;
+                           << "   ";
+                       
+                                                        // compute coordinates for
+                                                        // this patch point
+                       const Point<spacedim>
+                         node = ((((patch->vertices[1] * x_frac) +
+                                   (patch->vertices[0] * (1-x_frac))) * (1-y_frac) +
+                                  ((patch->vertices[2] * x_frac) +
+                                   (patch->vertices[3] * (1-x_frac))) * y_frac)   * (1-z_frac) +
+                                 (((patch->vertices[5] * x_frac) +
+                                   (patch->vertices[4] * (1-x_frac))) * (1-y_frac) +
+                                  ((patch->vertices[6] * x_frac) +
+                                   (patch->vertices[7] * (1-x_frac))) * y_frac)   * z_frac);
+                       
+                                                        // write out coordinates
+                       for (unsigned int i=0; i<spacedim; ++i)
+                         out << node(i) << ' ';
+                                                        // fill with zeroes unnecessary here
+                       for (unsigned int i=spacedim; i<3; ++i)
+                         out << "0 ";
+                       out << endl;
                        
                        ++present_node;
                      };
@@ -212,7 +248,7 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
       unsigned int present_cell = 1;      
       unsigned int first_vertex_of_patch = 0;
       
-      for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+      for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
           patch!=patches.end(); ++patch)
        {
          const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -276,7 +312,7 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
            };
 
 
-// finally update the number
+                                          // finally update the number
                                           // of the first vertex of this patch
          switch (dim)
            {
@@ -304,7 +340,7 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
     };
 
 
-/////////////////////////////
+                                  /////////////////////////////
                                   // now write data
   if (n_data_sets != 0)
     {      
@@ -320,9 +356,9 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
            << endl;
 
 
-// loop over all patches
+                                      // loop over all patches
       unsigned int present_node = 1;
-      for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+      for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
           patch != patches.end(); ++patch)
        {
          const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -407,11 +443,12 @@ void DataOutBase::write_ucd (const vector<Patch<dim> > &patches,
 };
 
 
-template <int dim>
-void DataOutBase::write_gnuplot (const vector<Patch<dim> > &patches,
-                                const vector<string>      &data_names,
-                                const GnuplotFlags        &/*flags*/,
-                                ostream                   &out) 
+
+template <int dim, int spacedim>
+void DataOutBase::write_gnuplot (const vector<Patch<dim,spacedim> > &patches,
+                                const vector<string>               &data_names,
+                                const GnuplotFlags                 &/*flags*/,
+                                ostream                            &out) 
 {
   AssertThrow (out, ExcIO());
   
@@ -442,7 +479,7 @@ void DataOutBase::write_gnuplot (const vector<Patch<dim> > &patches,
          << "#" << endl
          << "# ";
       
-      switch (dim) 
+      switch (spacedim)
        {
          case 1:
                out << "<x> ";
@@ -466,8 +503,8 @@ void DataOutBase::write_gnuplot (const vector<Patch<dim> > &patches,
     };
 
 
-// loop over all patches
-  for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+                                  // loop over all patches
+  for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
        patch != patches.end(); ++patch)
     {
       const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -565,7 +602,7 @@ void DataOutBase::write_gnuplot (const vector<Patch<dim> > &patches,
 
                                                     // compute coordinates for
                                                     // this patch point
-                   const Point<dim> this_point
+                   const Point<spacedim> this_point
                      = ((((patch->vertices[1] * x_frac) +
                           (patch->vertices[0] * (1-x_frac))) * (1-y_frac) +
                          ((patch->vertices[2] * x_frac) +
@@ -694,16 +731,18 @@ void DataOutBase::write_gnuplot (const vector<Patch<dim> > &patches,
 };
 
 
-template <int dim>
-void DataOutBase::write_povray (const vector<Patch<dim> > &patches,
-                               const vector<string>      &data_names,
-                               const PovrayFlags         &flags,
-                               ostream                   &out) 
+
+template <int dim, int spacedim>
+void DataOutBase::write_povray (const vector<Patch<dim,spacedim> > &patches,
+                               const vector<string>               &data_names,
+                               const PovrayFlags                  &flags,
+                               ostream                            &out) 
 {
   AssertThrow (out, ExcIO());
   
   Assert (patches.size() > 0, ExcNoPatches());
-  Assert (dim==2, ExcNotImplemented());        // only for 2-D  
+  Assert (dim==2, ExcNotImplemented());        // only for 2-D surfaces on a 2-D plane
+  Assert (spacedim==2, ExcNotImplemented());
 
   const unsigned int n_data_sets = data_names.size();
   
@@ -728,26 +767,26 @@ void DataOutBase::write_povray (const vector<Patch<dim> > &patches,
          << endl
          << "*/ " << endl;
       
-                                    // include files
+                                      // include files
       out << "#include \"colors.inc\" " << endl
          << "#include \"textures.inc\" " << endl;
 
-
-// use external include file for textures,
-                                    // camera and light
+      
+                                      // use external include file for textures,
+                                      // camera and light
       if (flags.external_data)
        out << "#include \"data.inc\" " << endl;
       else                          // all definitions in data file
        {  
-                                    // camera
+                                          // camera
          out << endl << endl
              << "camera {"            << endl
              << "  location <1,4,-7>" << endl
              << "  look_at <0,0,0>"   << endl
              << "  angle 30"          << endl
              << "}"                   << endl;
-      
-                                    // light
+         
+                                          // light
          out << endl 
              << "light_source {"      << endl
              << "  <1,4,-7>"      << endl
@@ -760,11 +799,11 @@ void DataOutBase::write_povray (const vector<Patch<dim> > &patches,
              << "}"                   << endl;
        }
     };
-
+  
                                   // max. and min. heigth of solution 
   double hmin=0, hmax=0;
 
-  for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+  for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
        patch != patches.end(); ++patch)
     {
       const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -814,7 +853,7 @@ void DataOutBase::write_povray (const vector<Patch<dim> > &patches,
     }
 
                                      // loop over all patches
-  for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+  for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
        patch != patches.end(); ++patch)
     {
       const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -831,7 +870,7 @@ void DataOutBase::write_povray (const vector<Patch<dim> > &patches,
              ExcInvalidDatasetSize (patch->data.n(), n_subdivisions+1));
 
 
-      Point<dim> ver[16];                          // value for all points in this patch
+      Point<spacedim> ver[16];                     // value for all points in this patch
       
       for (unsigned int i=0; i<n_subdivisions+1; ++i)
        {
@@ -994,15 +1033,17 @@ void DataOutBase::write_povray (const vector<Patch<dim> > &patches,
 };
 
 
-template <int dim>
-void DataOutBase::write_eps (const vector<Patch<dim> > &patches,
-                            const vector<string>      &/*data_names*/,
-                            const EpsFlags            &flags,
-                            ostream                   &out) 
+
+template <int dim, int spacedim>
+void DataOutBase::write_eps (const vector<Patch<dim,spacedim> > &patches,
+                            const vector<string>               &/*data_names*/,
+                            const EpsFlags                     &flags,
+                            ostream                            &out) 
 {
   Assert (out, ExcIO());
   
   Assert (patches.size() > 0, ExcNoPatches());
+  Assert (spacedim==dim, ExcNotImplemented());
   
   switch (dim) 
     {
@@ -1035,7 +1076,7 @@ void DataOutBase::write_eps (const vector<Patch<dim> > &patches,
                                         // note that since dim==2, we
                                         // have exactly four vertices per
                                         // patch and per cell
-       for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+       for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
             patch!=patches.end(); ++patch)
          {
            const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -1048,7 +1089,7 @@ void DataOutBase::write_eps (const vector<Patch<dim> > &patches,
                              x_frac1 = (i+1) * 1./n_subdivisions,
                              y_frac1 = (j+1) * 1./n_subdivisions;
                  
-                 const Point<dim> points[4]
+                 const Point<spacedim> points[4]
                    = { (((patch->vertices[1] * x_frac) +
                          (patch->vertices[0] * (1-x_frac))) * (1-y_frac) +
                         ((patch->vertices[2] * x_frac) +
@@ -1141,7 +1182,7 @@ void DataOutBase::write_eps (const vector<Patch<dim> > &patches,
 
                                                   // compute coordinates of
                                                   // center of cell
-                 const Point<dim> center_point
+                 const Point<spacedim> center_point
                    = (points[0] + points[1] + points[2] + points[3]) / 4;
                  const double center_height
                    = -(heights[0] + heights[1] + heights[2] + heights[3]) / 4;
@@ -1331,11 +1372,12 @@ void DataOutBase::write_eps (const vector<Patch<dim> > &patches,
 };
 
 
-template <int dim>
-void DataOutBase::write_gmv (const vector<Patch<dim> > &patches,
-                            const vector<string>      &data_names,
-                            const GmvFlags            &/*flags*/,
-                            ostream                   &out) 
+
+template <int dim, int spacedim>
+void DataOutBase::write_gmv (const vector<Patch<dim,spacedim> > &patches,
+                            const vector<string>               &data_names,
+                            const GmvFlags                     &/*flags*/,
+                            ostream                            &out) 
 {
   AssertThrow (out, ExcIO());
 
@@ -1360,7 +1402,7 @@ void DataOutBase::write_gmv (const vector<Patch<dim> > &patches,
                                   // and cells for later use
   unsigned int n_cells = 0,
               n_nodes = 0;
-  for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+  for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
        patch!=patches.end(); ++patch)
     switch (dim)
       {
@@ -1413,7 +1455,7 @@ void DataOutBase::write_gmv (const vector<Patch<dim> > &patches,
   Threads::ThreadManager thread_manager;
   Threads::spawn (thread_manager,
                  Threads::encapsulate (&DataOutBase::template
-                                       write_gmv_reorder_data_vectors<dim>)
+                                       write_gmv_reorder_data_vectors<dim,spacedim>)
                  .collect_args(patches, data_vectors));
 
                                   ///////////////////////////////
@@ -1426,14 +1468,14 @@ void DataOutBase::write_gmv (const vector<Patch<dim> > &patches,
   out << "nodes " << n_nodes << endl;
   for (unsigned int d=1; d<=3; ++d)
     {
-      for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+      for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
           patch!=patches.end(); ++patch)
        {
          const unsigned int n_subdivisions = patch->n_subdivisions;
          
                                           // if we have nonzero values for
                                           // this coordinate
-         if (d<=dim)
+         if (d<=spacedim)
            {
              switch (dim)
                {
@@ -1503,7 +1545,7 @@ void DataOutBase::write_gmv (const vector<Patch<dim> > &patches,
                };
            }
          else
-                                            // d>dim. write zeros instead
+                                            // d>spacedim. write zeros instead
            {
              const unsigned int n_points
                = static_cast<unsigned int>(pow (n_subdivisions+1, dim));
@@ -1526,7 +1568,7 @@ void DataOutBase::write_gmv (const vector<Patch<dim> > &patches,
 
       unsigned int first_vertex_of_patch = 0;
       
-      for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+      for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
           patch!=patches.end(); ++patch)
        {
          const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -1644,9 +1686,10 @@ void DataOutBase::write_gmv (const vector<Patch<dim> > &patches,
 
 
 
-template <int dim>
-void DataOutBase::write_gmv_reorder_data_vectors (const vector<Patch<dim> > &patches,
-                                                 vector<vector<double> >   &data_vectors)
+template <int dim, int spacedim>
+void
+DataOutBase::write_gmv_reorder_data_vectors (const vector<Patch<dim,spacedim> > &patches,
+                                            vector<vector<double> >            &data_vectors)
 {
                                   // unlike in the main function, we
                                   // don't have here the data_names
@@ -1662,7 +1705,7 @@ void DataOutBase::write_gmv_reorder_data_vectors (const vector<Patch<dim> > &pat
   
                                   // loop over all patches
   unsigned int next_value = 0;
-  for (typename vector<Patch<dim> >::const_iterator patch=patches.begin();
+  for (typename vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
        patch != patches.end(); ++patch)
     {
       const unsigned int n_subdivisions = patch->n_subdivisions;
@@ -1735,48 +1778,53 @@ void DataOutBase::write_gmv_reorder_data_vectors (const vector<Patch<dim> > &pat
 /* --------------------------- class DataOutInterface ---------------------- */
 
 
-template <int dim>
-void DataOutInterface<dim>::write_ucd (ostream &out) const 
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write_ucd (ostream &out) const 
 {
   DataOutBase::write_ucd (get_patches(), get_dataset_names(),
                          ucd_flags, out);
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::write_gnuplot (ostream &out) const 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write_gnuplot (ostream &out) const 
 {
   DataOutBase::write_gnuplot (get_patches(), get_dataset_names(),
                              gnuplot_flags, out);
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::write_povray (ostream &out) const 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write_povray (ostream &out) const 
 {
   DataOutBase::write_povray (get_patches(), get_dataset_names(),
                             povray_flags, out);
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::write_eps (ostream &out) const 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write_eps (ostream &out) const 
 {
   DataOutBase::write_eps (get_patches(), get_dataset_names(),
                          eps_flags, out);
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::write_gmv (ostream &out) const 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write_gmv (ostream &out) const 
 {
   DataOutBase::write_gmv (get_patches(), get_dataset_names(),
                          gmv_flags, out);
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::write (ostream &out,
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write (ostream &out,
                                   OutputFormat output_format) const
 {
   if (output_format == default_format)
@@ -1810,51 +1858,58 @@ void DataOutInterface<dim>::write (ostream &out,
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::set_default_format(OutputFormat fmt)
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_default_format(OutputFormat fmt)
 {
   Assert(fmt != default_format, ExcNotImplemented());
   default_fmt = fmt;
 }
 
 
-template <int dim>
-void DataOutInterface<dim>::set_flags (const UcdFlags &flags) 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_flags (const UcdFlags &flags) 
 {
   ucd_flags = flags;
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::set_flags (const GnuplotFlags &flags) 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_flags (const GnuplotFlags &flags) 
 {
   gnuplot_flags = flags;
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::set_flags (const PovrayFlags &flags) 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_flags (const PovrayFlags &flags) 
 {
   povray_flags = flags;
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::set_flags (const EpsFlags &flags) 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_flags (const EpsFlags &flags) 
 {
   eps_flags = flags;
 };
 
 
-template <int dim>
-void DataOutInterface<dim>::set_flags (const GmvFlags &flags) 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_flags (const GmvFlags &flags) 
 {
   gmv_flags = flags;
 };
 
 
-template <int dim>
-string DataOutInterface<dim>::default_suffix (OutputFormat output_format) const
+
+template <int dim, int spacedim>
+string DataOutInterface<dim,spacedim>::default_suffix (OutputFormat output_format) const
 {
   if (output_format == default_format)
     output_format = default_fmt;
@@ -1883,9 +1938,10 @@ string DataOutInterface<dim>::default_suffix (OutputFormat output_format) const
 };
 
 
-template <int dim>
-DataOutInterface<dim>::OutputFormat
-DataOutInterface<dim>::parse_output_format (const string &format_name)
+
+template <int dim, int spacedim>
+DataOutInterface<dim,spacedim>::OutputFormat
+DataOutInterface<dim,spacedim>::parse_output_format (const string &format_name)
 {
   if (format_name == "ucd")
     return ucd;
@@ -1909,15 +1965,17 @@ DataOutInterface<dim>::parse_output_format (const string &format_name)
 };
 
 
-template <int dim>
-string DataOutInterface<dim>::get_output_format_names ()
+
+template <int dim, int spacedim>
+string DataOutInterface<dim,spacedim>::get_output_format_names ()
 {
   return "ucd|gnuplot|povray|eps|gmv";
 }
 
 
-template <int dim>
-void DataOutInterface<dim>::declare_parameters (ParameterHandler &prm) 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::declare_parameters (ParameterHandler &prm) 
 {
   prm.declare_entry ("Output format", "gnuplot",
                     Patterns::Selection (get_output_format_names ()));
@@ -1944,8 +2002,9 @@ void DataOutInterface<dim>::declare_parameters (ParameterHandler &prm)
 }
 
 
-template <int dim>
-void DataOutInterface<dim>::parse_parameters (ParameterHandler &prm) 
+
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::parse_parameters (ParameterHandler &prm) 
 {
   const string& output_name = prm.get ("Output format");
   default_fmt = parse_output_format (output_name);
@@ -1974,5 +2033,11 @@ void DataOutInterface<dim>::parse_parameters (ParameterHandler &prm)
 
 // explicit instantiations. functions in DataOutBase are instantiated by
 // the respective functions in DataOut_Interface
-template class DataOutInterface<data_out_dimension>;
-template class DataOutBase::Patch<data_out_dimension>;
+template class DataOutInterface<data_out_dimension,data_out_dimension>;
+template class DataOutBase::Patch<data_out_dimension, data_out_dimension>;
+
+// also enable plotting surfaces of 3d objects
+#if data_out_dimension == 3
+template class DataOutInterface<data_out_dimension,data_out_dimension-1>;
+template class DataOutBase::Patch<data_out_dimension, data_out_dimension-1>;
+#endif
index 37f389869f251f10ab911f4cadaf3faf065a6038..18ab0cbad7aa884f8d9ff5fbe90555a70a5b5186 100644 (file)
@@ -9,8 +9,22 @@
     <meta name="keywords" content="deal.II"></head>
 
 <body>
+
+
 <h2>Changes after Version 3.0</h2>
+
+This is a quite extensive list of changes made after the release of 
+<acronym>deal.II</acronym> version 3.0. It is subdivided into changes
+made to the three sub-libraries <a href="#base">base</a>, 
+<a href="#lac">lac</a>, and <a href="#deal.II">deal.II</a>, as well as
+changes to the <a href="#general">general infrastructure,
+documentation, etc</a>.
+
+
+
+<a name="general"></a>
 <h3>General</h3>
+
 <ol>
   <li> <p>
        New: the <code class="program">step-9</code> example program is
        <br>
        (GK 2000/04/05)
        </p>
-
 </ol>
 
+
+
+<a name="base"></a>
 <h3>base</h3>
 
 <ol>
+  <li> <p>
+       New: The <code class="class">DataOutBase</code> and
+       <code class="class">DataOutBase::Patch</code> classes have been
+       changed so as to allow output of objects that have an other
+       dimension than the surrounding space, for example writing faces
+       instead of cells (this might be useful to write only external
+       faces in 3d computations).
+       <br>
+       (WB 2000/09/07)
+       </p>       
+
   <li> <p> 
        New: <code class="class">Timer</code> now uses the system
        function <code class="member">getrusage (RUSAGE_CHILDREN,
 
 
 
+<a name="lac"></a>
 <h3>lac</h3>
 
 <ol>
 
 
 
+<a name="deal.II"></a>
 <h3>deal.II</h3>
 
 <ol>

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.