* details of the viewpoint, light source, etc.
*
*
- * @author Wolfgang Bangerth 1999, 2000, 2001; 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, 2001; EPS output based on an earlier implementation by Stefan Nauber for the old DataOut class; Povray output by Thomas Richter 1999; OpenDX output by Guido Kanschat, 2001; Tecplot output by Benjamin Shelton Kirk, 2002.
*/
class DataOutBase
{
* present no flags are implemented.
*/
struct GmvFlags
+ {
+ private:
+ /**
+ * Dummy entry to suppress compiler
+ * warnings when copying an empty
+ * structure. Remove this member
+ * when adding the first flag to
+ * this structure (and remove the
+ * @p{private} as well).
+ */
+ int dummy;
+
+ public:
+ /**
+ * Declare all flags with name
+ * and type as offered by this
+ * class, for use in input files.
+ */
+ static void declare_parameters (ParameterHandler &prm);
+
+ /**
+ * Read the parameters declared in
+ * @p{declare_parameters} and set the
+ * flags for this output format
+ * accordingly.
+ *
+ * The flags thus obtained overwrite
+ * all previous contents of this object.
+ */
+ void parse_parameters (ParameterHandler &prm);
+
+ /**
+ * Determine an estimate for
+ * the memory consumption (in
+ * bytes) of this
+ * object. Since sometimes
+ * the size of objects can
+ * not be determined exactly
+ * (for example: what is the
+ * memory consumption of an
+ * STL @p{std::map} type with a
+ * certain number of
+ * elements?), this is only
+ * an estimate. however often
+ * quite close to the true
+ * value.
+ */
+ unsigned int memory_consumption () const;
+ };
+
+ /**
+ * Flags describing the details of
+ * output in Tecplot format. At
+ * present no flags are implemented.
+ */
+ struct TecplotFlags
{
private:
/**
const GmvFlags &flags,
std::ostream &out);
+ /**
+ * Write the given list of patches
+ * to the output stream in Tecplot
+ * format. See the general
+ * documentation for more information
+ * on the parameters.
+ */
+ template <int dim, int spacedim>
+ static void write_tecplot (const typename std::vector<Patch<dim,spacedim> > &patches,
+ const std::vector<std::string> &data_names,
+ const TecplotFlags &flags,
+ std::ostream &out);
+
/**
* Write the given list of
* patches to the output stream
* the presently supported output
* formats.
*/
- enum OutputFormat { default_format, dx, ucd, gnuplot, povray, eps, gmv, vtk };
+ enum OutputFormat { default_format, dx, ucd, gnuplot, povray, eps, gmv, tecplot, vtk };
/**
* Obtain data through the
*/
void write_gmv (std::ostream &out) const;
+
+ /**
+ * Obtain data through the
+ * @p{get_patches} function and
+ * write it to @p{out} in Tecplot
+ * format.
+ */
+ void write_tecplot (std::ostream &out) const;
+
/**
* Obtain data through the
* @p{get_patches} function and
*/
void set_flags (const GmvFlags &gmv_flags);
+ /**
+ * Set the flags to be used for
+ * output in Tecplot format.
+ */
+ void set_flags (const TecplotFlags &tecplot_flags);
+
/**
* Set the flags to be used for
* output in GMV format.
*/
GmvFlags gmv_flags;
+ /**
+ * Flags to be used upon output
+ * of Tecplot data in one space
+ * dimension. Can be changed by
+ * using the @p{set_flags}
+ * function.
+ */
+ TecplotFlags tecplot_flags;
+
/**
* Flags to be used upon output
* of vtk data in one space
+template <int dim, int spacedim>
+void DataOutBase::write_tecplot (const typename std::vector<Patch<dim,spacedim> > &patches,
+ const std::vector<std::string> &data_names,
+ const TecplotFlags &/*flags*/,
+ std::ostream &out)
+{
+ AssertThrow (out, ExcIO());
+
+ Assert (patches.size() > 0, ExcNoPatches());
+
+ const unsigned int n_data_sets = data_names.size();
+ // check against # of data sets in
+ // first patch. checks against all
+ // other patches are made in
+ // write_gmv_reorder_data_vectors
+ Assert (n_data_sets == patches[0].data.n_rows(),
+ ExcUnexpectedNumberOfDatasets (patches[0].data.n_rows(), n_data_sets));
+
+
+
+ // first count the number of cells
+ // and cells for later use
+ unsigned int n_cells = 0,
+ n_nodes = 0;
+ for (typename std::vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
+ patch!=patches.end(); ++patch)
+ switch (dim)
+ {
+ case 1:
+ n_cells += patch->n_subdivisions;
+ n_nodes += patch->n_subdivisions+1;
+ break;
+ case 2:
+ n_cells += patch->n_subdivisions *
+ patch->n_subdivisions;
+ n_nodes += (patch->n_subdivisions+1) *
+ (patch->n_subdivisions+1);
+ break;
+ case 3:
+ n_cells += patch->n_subdivisions *
+ patch->n_subdivisions *
+ patch->n_subdivisions;
+ n_nodes += (patch->n_subdivisions+1) *
+ (patch->n_subdivisions+1) *
+ (patch->n_subdivisions+1);
+ break;
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+
+
+ ///////////////////////
+ // preamble
+ {
+ std::time_t time1= std::time (0);
+ std::tm *time = std::localtime(&time1);
+ out << "# This file was generated by the deal.II library." << std::endl
+ << "# Date = "
+ << time->tm_year+1900 << "/"
+ << time->tm_mon+1 << "/"
+ << time->tm_mday << std::endl
+ << "# Time = "
+ << time->tm_hour << ":"
+ << std::setw(2) << time->tm_min << ":"
+ << std::setw(2) << time->tm_sec << std::endl
+ << "#" << std::endl
+ << "# For a description of the Tecplot format see the Tecplot documentation."
+ << std::endl
+ << "#" << std::endl;
+
+
+ out << "Variables=";
+
+ switch (dim)
+ {
+ case 1:
+ out << "\"x\"";
+ break;
+ case 2:
+ out << "\"x\", \"y\"";
+ break;
+ case 3:
+ out << "\"x\", \"y\", \"z\"";
+ break;
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+
+ for (unsigned int data_set=0; data_set<n_data_sets; ++data_set)
+ {
+ out << ", \"" << data_names[data_set] << "\"";
+ };
+
+ out << std::endl;
+
+ if (dim > 1)
+ {
+ out << "zone f=feblock, n=" << n_nodes << ", e=" << n_cells << ", et=";
+
+ switch (dim)
+ {
+ case 2:
+ out << "quadrilateral" << std::endl;
+ break;
+ case 3:
+ out << "brick" << std::endl;
+ break;
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+ }
+ else
+ {
+ out << "zone f=block, n=" << n_nodes << std::endl;
+ };
+ };
+
+
+ // in Tecplot block format the vertex
+ // coordinates and the data have an
+ // order that is a bit unpleasant
+ // (first all x coordinates, then
+ // all y coordinate, ...; first all
+ // data of variable 1, then
+ // variable 2, etc), so we have to
+ // copy the data vectors a bit around
+ //
+ // note that we copy vectors when
+ // looping over the patches since we
+ // have to write them one variable
+ // at a time and don't want to use
+ // more than one loop
+ //
+ // this copying of data vectors can
+ // be done while we already output
+ // the vertices, so do this on a
+ // separate thread and when wanting
+ // to write out the data, we wait
+ // for that thread to finish
+ std::vector<std::vector<double> > data_vectors (n_data_sets,
+ std::vector<double> (n_nodes));
+ Threads::ThreadManager thread_manager;
+ void (*fun_ptr) (const typename std::vector<Patch<dim,spacedim> > &,
+ std::vector<std::vector<double> > &)
+ = &DataOutBase::template write_gmv_reorder_data_vectors<dim,spacedim>;
+ Threads::spawn (thread_manager,
+ Threads::encapsulate (fun_ptr)
+ .collect_args(patches, data_vectors));
+
+ ///////////////////////////////
+ // first make up a list of used
+ // vertices along with their
+ // coordinates
+
+
+ for (unsigned int d=1; d<=spacedim; ++d)
+ {
+ for (typename std::vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
+ patch!=patches.end(); ++patch)
+ {
+ const unsigned int n_subdivisions = patch->n_subdivisions;
+
+ switch (dim)
+ {
+ case 1:
+ {
+ for (unsigned int i=0; i<n_subdivisions+1; ++i)
+ out << ((patch->vertices[1](0) * i / n_subdivisions) +
+ (patch->vertices[0](0) * (n_subdivisions-i) / n_subdivisions))
+ << '\n';
+ break;
+ };
+
+ case 2:
+ {
+ for (unsigned int i=0; i<n_subdivisions+1; ++i)
+ for (unsigned int j=0; j<n_subdivisions+1; ++j)
+ {
+ const double x_frac = i * 1./n_subdivisions,
+ y_frac = j * 1./n_subdivisions;
+
+ // compute coordinates for
+ // this patch point
+ out << (((patch->vertices[1](d-1) * x_frac) +
+ (patch->vertices[0](d-1) * (1-x_frac))) * (1-y_frac) +
+ ((patch->vertices[2](d-1) * x_frac) +
+ (patch->vertices[3](d-1) * (1-x_frac))) * y_frac)
+ << '\n';
+ };
+ break;
+ };
+
+ case 3:
+ {
+ for (unsigned int i=0; i<n_subdivisions+1; ++i)
+ for (unsigned int j=0; j<n_subdivisions+1; ++j)
+ for (unsigned int k=0; k<n_subdivisions+1; ++k)
+ {
+ // note the broken
+ // design of hexahedra
+ // in deal.II, where
+ // first the z-component
+ // is counted up, before
+ // increasing the y-
+ // coordinate
+ const double x_frac = i * 1./n_subdivisions,
+ y_frac = k * 1./n_subdivisions,
+ z_frac = j * 1./n_subdivisions;
+
+ // compute coordinates for
+ // this patch point
+ out << ((((patch->vertices[1](d-1) * x_frac) +
+ (patch->vertices[0](d-1) * (1-x_frac))) * (1-y_frac) +
+ ((patch->vertices[2](d-1) * x_frac) +
+ (patch->vertices[3](d-1) * (1-x_frac))) * y_frac) * (1-z_frac) +
+ (((patch->vertices[5](d-1) * x_frac) +
+ (patch->vertices[4](d-1) * (1-x_frac))) * (1-y_frac) +
+ ((patch->vertices[6](d-1) * x_frac) +
+ (patch->vertices[7](d-1) * (1-x_frac))) * y_frac) * z_frac)
+ << '\n';
+ };
+
+ break;
+ };
+
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+ };
+ out << std::endl;
+ };
+
+
+ ///////////////////////////////////////
+ // data output.
+ //
+ // now write the data vectors to
+ // @p{out} first make sure that all
+ // data is in place
+ thread_manager.wait ();
+
+ // then write data.
+ for (unsigned int data_set=0; data_set<n_data_sets; ++data_set)
+ {
+ copy(data_vectors[data_set].begin(),
+ data_vectors[data_set].end(),
+ std::ostream_iterator<double>(out, "\n"));
+ out << std::endl;
+ };
+
+
+ /////////////////////////////////
+ // now for the cells. note that
+ // vertices are counted from 1 onwards
+
+ unsigned int first_vertex_of_patch = 0;
+
+ for (typename std::vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
+ patch!=patches.end(); ++patch)
+ {
+ const unsigned int n_subdivisions = patch->n_subdivisions;
+
+ // write out the cells making
+ // up this patch
+ switch (dim)
+ {
+ case 1:
+ {
+ break;
+ };
+ case 2:
+ {
+ for (unsigned int i=0; i<n_subdivisions; ++i)
+ for (unsigned int j=0; j<n_subdivisions; ++j)
+ out << first_vertex_of_patch+i*(n_subdivisions+1)+j+1 << ' '
+ << first_vertex_of_patch+(i+1)*(n_subdivisions+1)+j+1 << ' '
+ << first_vertex_of_patch+(i+1)*(n_subdivisions+1)+j+1+1 << ' '
+ << first_vertex_of_patch+i*(n_subdivisions+1)+j+1+1
+ << std::endl;
+ break;
+ };
+
+ case 3:
+ {
+ for (unsigned int i=0; i<n_subdivisions; ++i)
+ for (unsigned int j=0; j<n_subdivisions; ++j)
+ for (unsigned int k=0; k<n_subdivisions; ++k)
+ {
+ // note: vertex indices start with 1!
+ out << first_vertex_of_patch+(i*(n_subdivisions+1)+j )*(n_subdivisions+1)+k +1 << ' '
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j )*(n_subdivisions+1)+k +1 << ' '
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j+1)*(n_subdivisions+1)+k +1 << ' '
+ << first_vertex_of_patch+(i*(n_subdivisions+1)+j+1 )*(n_subdivisions+1)+k +1 << ' '
+ << first_vertex_of_patch+(i*(n_subdivisions+1)+j )*(n_subdivisions+1)+k+1+1 << ' '
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j )*(n_subdivisions+1)+k+1+1 << ' '
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j+1)*(n_subdivisions+1)+k+1+1 << ' '
+ << first_vertex_of_patch+(i*(n_subdivisions+1)+j+1 )*(n_subdivisions+1)+k+1+1 << ' '
+ << std::endl;
+ };
+ break;
+ };
+
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+
+
+ // finally update the number
+ // of the first vertex of this patch
+ switch (dim)
+ {
+ case 1:
+ first_vertex_of_patch += n_subdivisions+1;
+ break;
+ case 2:
+ first_vertex_of_patch += (n_subdivisions+1) *
+ (n_subdivisions+1);
+ break;
+ case 3:
+ first_vertex_of_patch += (n_subdivisions+1) *
+ (n_subdivisions+1) *
+ (n_subdivisions+1);
+ break;
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+ };
+
+ // assert the stream is still ok
+ AssertThrow (out, ExcIO());
+};
+
+
+
template <int dim, int spacedim>
void DataOutBase::write_vtk (const typename std::vector<Patch<dim,spacedim> > &patches,
const std::vector<std::string> &data_names,
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write_tecplot (std::ostream &out) const
+{
+ DataOutBase::write_tecplot (get_patches(), get_dataset_names(),
+ tecplot_flags, out);
+};
+
+
+
template <int dim, int spacedim>
void DataOutInterface<dim,spacedim>::write_vtk (std::ostream &out) const
{
write_gmv (out);
break;
+ case tecplot:
+ write_tecplot (out);
+ break;
+
case vtk:
write_vtk (out);
break;
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_flags (const TecplotFlags &flags)
+{
+ tecplot_flags = flags;
+};
+
+
+
template <int dim, int spacedim>
void DataOutInterface<dim,spacedim>::set_flags (const VtkFlags &flags)
{
case gmv:
return ".gmv";
+
+ case tecplot:
+ return ".dat";
case vtk:
return ".vtk";
if (format_name == "gmv")
return gmv;
+
+ if (format_name == "tecplot")
+ return tecplot;
if (format_name == "vtk")
return vtk;
template <int dim, int spacedim>
std::string DataOutInterface<dim,spacedim>::get_output_format_names ()
{
- return "dx|ucd|gnuplot|povray|eps|gmv|vtk";
+ return "dx|ucd|gnuplot|povray|eps|gmv|tecplot|vtk";
}
GmvFlags::declare_parameters (prm);
prm.leave_subsection ();
+ prm.enter_subsection ("Tecplot output parameters");
+ TecplotFlags::declare_parameters (prm);
+ prm.leave_subsection ();
+
prm.enter_subsection ("Vtk output parameters");
VtkFlags::declare_parameters (prm);
prm.leave_subsection ();
gmv_flags.parse_parameters (prm);
prm.leave_subsection ();
+ prm.enter_subsection ("Tecplot output parameters");
+ tecplot_flags.parse_parameters (prm);
+ prm.leave_subsection ();
+
prm.enter_subsection ("Vtk output parameters");
vtk_flags.parse_parameters (prm);
prm.leave_subsection ();
MemoryConsumption::memory_consumption (povray_flags) +
MemoryConsumption::memory_consumption (eps_flags) +
MemoryConsumption::memory_consumption (gmv_flags) +
+ MemoryConsumption::memory_consumption (tecplot_flags) +
MemoryConsumption::memory_consumption (vtk_flags));
};