]> https://gitweb.dealii.org/ - dealii.git/commitdiff
enable filename selection
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 22 May 2014 22:05:31 +0000 (22:05 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 22 May 2014 22:05:31 +0000 (22:05 +0000)
git-svn-id: https://svn.dealii.org/trunk@32970 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/numerics/dof_output_operator.h
deal.II/include/deal.II/numerics/dof_output_operator.templates.h

index 70f89b546791d6cbad742b0e567c4c7b14d06494..b9f55a04a21ccdc043e6abc42e42de6e8a8c4751 100644 (file)
@@ -30,18 +30,37 @@ DEAL_II_NAMESPACE_OPEN
 
 namespace Algorithms
 {
+  /**
+   * An output operator writing a separate file in each step and
+   * writing the vectors as finite element functions with respect to a
+   * given DoFHandler.
+   */
   template <class VECTOR, int dim, int spacedim=dim>
   class DoFOutputOperator : public OutputOperator<VECTOR>
   {
   public:
+      /*
+       * Constructor. The <tt>filename</tt> is the common base name of all
+       * files and the argument <tt>digits</tt> should be the number of digits
+       * of the highest number in the sequence. File names by default
+       * have the form "outputNN" with NNN the number set by the last
+       * step command. Numbers with less digits are filled with zeros
+       * from the left.
+       */
+      DoFOutputOperator (const std::string filename_base = std::string("output"),
+      const unsigned int digits = 3);
+      
     void initialize (DoFHandler<dim, spacedim> &dof_handler);
 
     virtual OutputOperator<VECTOR> &
-    operator<<(const NamedData<VECTOR *> &vectors);
+    operator << (const AnyData &vectors);
 
   private:
     SmartPointer<DoFHandler<dim, spacedim>,
                  DoFOutputOperator<VECTOR, dim, spacedim> > dof;
+
+      const std::string filename_base;
+      const unsigned int digits;
   };
 
   template <class VECTOR, int dim, int spacedim>
index b9529a01c48e0c1ec6bd14e8338859c90ef41d06..2bdb0f72e3d9a9dcd2af654e7c3160c63aa1be5e 100644 (file)
@@ -22,19 +22,37 @@ DEAL_II_NAMESPACE_OPEN
 
 namespace Algorithms
 {
+  template <class VECTOR, int dim, int spacedim>
+  DoFOutputOperator<VECTOR, dim, spacedim>::DoFOutputOperator (
+    const std::string filename_base,
+    const unsigned int digits)
+                 :
+                 filename_base(filename_base),
+                 digits(digits)
+  {}
+  
+
   template <class VECTOR, int dim, int spacedim>
   OutputOperator<VECTOR> &
   DoFOutputOperator<VECTOR, dim, spacedim>::operator<<(
-    const NamedData<VECTOR *> &vectors)
+    const AnyData &data)
   {
     Assert ((dof!=0), ExcNotInitialized());
     DataOut<dim> out;
+    out.set_default_format(DataOutBase::gnuplot);
     out.attach_dof_handler (*dof);
-    out.add_data_vector (*vectors(vectors.find("solution")), "solution");
-    out.add_data_vector (*vectors(vectors.find("update")), "update");
-    out.add_data_vector (*vectors(vectors.find("residual")), "residual");
+    for (unsigned int i=0;i<data.size();++i)
+      {
+       const VECTOR* p = data.try_read_ptr<VECTOR>(i);
+       if (p!=0)
+         {
+           out.add_data_vector (*p, data.name(i));
+         }
+      }
     std::ostringstream streamOut;
-    streamOut << "Newton_" << std::setw(3) << std::setfill('0') << this->step;
+    streamOut << filename_base
+             << std::setw(digits) << std::setfill('0') << this->step
+             << out.default_suffix();
     std::ofstream out_filename (streamOut.str().c_str());
     out.build_patches (2);
     out.write_gnuplot (out_filename);

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.