*
* @sect2{Supported output formats}
*
+ * @sect3{OpenDX (IBM Open Visualization Data Explorer}
+ *
+ * Since Data Explorer (DX) is distributed as OpenSource, there is a
+ * well-accessible visualization tool for all (at least Unix-based)
+ * platforms. Therefore, output in its natural file format is
+ * included. Right now, this output is very basic, since the
+ * functionality of @p{DataOutBase} and derived classes has to be
+ * enhanced to use all of its features. Still, using some
+ * sophisticated networks, most can be recovered.
+ *
+ *
* @sect3{UCD format}
*
* The UCD format is described in the AVS developer's guide. Due to
};
+ /**
+ * Flags describing the details of
+ * output in OpenDX format. At
+ * present no flags are implemented.
+ */
+ struct DXFlags
+ {
+ 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.
+ */
+ unsigned int memory_consumption () const;
+ };
+
/**
* Flags describing the details of
* output in UCD format.
unsigned int memory_consumption () const;
};
+ /**
+ * Write the given list of patches
+ * to the output stream in OpenDX
+ * format. See the general
+ * documentation for more information
+ * on the parameters.
+ */
+ template <int dim, int spacedim>
+ static void write_dx (const typename std::vector<Patch<dim,spacedim> > &patches,
+ const std::vector<std::string> &data_names,
+ const DXFlags &flags,
+ std::ostream &out);
+
/**
* Write the given list of patches
* to the output stream in ucd
* the presently supported output
* formats.
*/
- enum OutputFormat { default_format, ucd, gnuplot, povray, eps, gmv, vtk };
+ enum OutputFormat { default_format, dx, ucd, gnuplot, povray, eps, gmv, vtk };
+
+ /**
+ * Obtain data through the
+ * @p{get_patches} function and
+ * write it to @p{out} in OpenDX
+ * format.
+ */
+ void write_dx (std::ostream &out) const;
/**
* Obtain data through the
*/
void set_default_format (const OutputFormat default_format);
+ /**
+ * Set the flags to be used for
+ * output in OpenDX format.
+ */
+ void set_flags (const DXFlags &dx_flags);
+
/**
* Set the flags to be used for
* output in UCD format.
* usually has. At present the following
* formats are defined:
* @begin{itemize}
+ * @item @p{dx}: @p{.dx}
* @item @p{ucd}: @p{.inp}
* @item @p{gnuplot}: @p{.gnuplot}
* @item @p{povray}: @p{.pov}
*/
OutputFormat default_fmt;
+ /**
+ * Flags to be used upon output
+ * of OpenDX data. Can be changed by
+ * using the @p{set_flags}
+ * function.
+ */
+ DXFlags dx_flags;
+
/**
* Flags to be used upon output
* of UCD data. Can be changed by
+template <int dim, int spacedim>
+void DataOutBase::write_dx (const typename std::vector<Patch<dim,spacedim> > &patches,
+ const std::vector<std::string> &data_names,
+ const DXFlags &flags,
+ std::ostream &out)
+{
+ AssertThrow (out, ExcIO());
+
+ Assert (patches.size() > 0, ExcNoPatches());
+
+ const unsigned int n_data_sets = data_names.size();
+
+ // 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
+//TODO:[GK] Fill in preamble later
+ if (false)
+ {
+ 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 UCD format see the AVS Developer's guide."
+ << std::endl
+ << "#" << std::endl;
+ };
+
+ // start with vertices
+ out << "object 1 class array type float rank 1 shape " << spacedim << " items " << n_nodes << " data follows"
+ << std::endl;
+
+ ///////////////////////////////
+ // first write the coordinates of all vertices
+ if (true)
+ {
+ for (typename std::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
+ switch (dim)
+ {
+ case 1:
+ {
+ for (unsigned int i=0; i<n_subdivisions+1; ++i)
+ {
+ 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) << '\t';
+ out << std::endl;
+ };
+
+ 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
+ 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) << '\t';
+ out << std::endl;
+ };
+
+ 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
+ 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) << '\t';
+ out << std::endl;
+ };
+
+ break;
+ };
+
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+ };
+ };
+
+ /////////////////////////////////////////
+ // write cells
+ out << "object 2 class array type int rank 1 shape " << GeometryInfo<dim>::vertices_per_cell
+ << " items " << n_cells << " data follows"
+ << std::endl;
+
+ if (true)
+ {
+ 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:
+ {
+ for (unsigned int i=0; i<n_subdivisions; ++i)
+ out << first_vertex_of_patch+i << '\t'
+ << first_vertex_of_patch+i+1 << std::endl;
+ 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 << '\t'
+ << first_vertex_of_patch+(i+1)*(n_subdivisions+1)+j << '\t'
+ << first_vertex_of_patch+i*(n_subdivisions+1)+j+1 << '\t'
+ << first_vertex_of_patch+(i+1)*(n_subdivisions+1)+j+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)
+ {
+//TODO:[GK] Put in correct order
+ out << first_vertex_of_patch+(i*(n_subdivisions+1)+j )*(n_subdivisions+1)+k << '\t'
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j )*(n_subdivisions+1)+k << '\t'
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j+1)*(n_subdivisions+1)+k << '\t'
+ << first_vertex_of_patch+(i*(n_subdivisions+1)+j+1 )*(n_subdivisions+1)+k << '\t'
+ << first_vertex_of_patch+(i*(n_subdivisions+1)+j )*(n_subdivisions+1)+k+1 << '\t'
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j )*(n_subdivisions+1)+k+1 << '\t'
+ << first_vertex_of_patch+((i+1)*(n_subdivisions+1)+j+1)*(n_subdivisions+1)+k+1 << '\t'
+ << first_vertex_of_patch+(i*(n_subdivisions+1)+j+1 )*(n_subdivisions+1)+k+1 << '\t'
+ << 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());
+ };
+ };
+ out << std::endl;
+ };
+
+ out << "attribute \"element type\" string \"quads\"" << std::endl
+ << "attribute \"ref\" string \"positions\"" << std::endl;
+
+ /////////////////////////////
+ // now write data
+ if (n_data_sets != 0)
+ {
+ out << "object 3 class array type float rank 1 shape " << n_data_sets
+ << " items " << n_nodes << " data follows" << std::endl;
+
+ // loop over all patches
+ for (typename std::vector<Patch<dim,spacedim> >::const_iterator patch=patches.begin();
+ patch != patches.end(); ++patch)
+ {
+ const unsigned int n_subdivisions = patch->n_subdivisions;
+
+ Assert (patch->data.n_rows() == n_data_sets,
+ ExcUnexpectedNumberOfDatasets (patch->data.n_rows(), n_data_sets));
+ Assert (patch->data.n_cols() == (dim==1 ?
+ n_subdivisions+1 :
+ (dim==2 ?
+ (n_subdivisions+1)*(n_subdivisions+1) :
+ (dim==3 ?
+ (n_subdivisions+1)*(n_subdivisions+1)*(n_subdivisions+1) :
+ 0))),
+ ExcInvalidDatasetSize (patch->data.n_cols(), n_subdivisions+1));
+
+ switch (dim)
+ {
+ case 1:
+ {
+ for (unsigned int i=0; i<n_subdivisions+1; ++i)
+ {
+ for (unsigned int data_set=0; data_set<n_data_sets; ++data_set)
+ out << patch->data(data_set,i) << '\t';
+ out << std::endl;
+ };
+
+ break;
+ };
+
+ case 2:
+ {
+ for (unsigned int i=0; i<n_subdivisions+1; ++i)
+ for (unsigned int j=0; j<n_subdivisions+1; ++j)
+ {
+ for (unsigned int data_set=0; data_set<n_data_sets; ++data_set)
+ out << patch->data(data_set,i*(n_subdivisions+1) + j) << '\t';
+ out << std::endl;
+ };
+
+ 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)
+ {
+ for (unsigned int data_set=0; data_set<n_data_sets; ++data_set)
+ out << patch->data(data_set,
+ (i*(n_subdivisions+1)+j)*(n_subdivisions+1)+k)
+ << '\t';
+ out << std::endl;
+ };
+
+ break;
+ };
+
+ default:
+ Assert (false, ExcNotImplemented());
+ };
+ };
+ out << "attribute \"dep\" string \"positions\"" << std::endl;
+ } else {
+ out << "object 3 class constantarray type float rank 0 items " << n_nodes << " data follows"
+ << std::endl << '0' << std::endl;
+ }
+
+ // no model data
+
+ out << "object \"deal data\" class field" << std::endl
+ << "component \"positions\" value 1" << std::endl
+ << "component \"connections\" value 2" << std::endl
+ << "component \"data\" value 3" << std::endl;
+
+ out << "end" << std::endl;
+ // assert the stream is still ok
+ AssertThrow (out, ExcIO());
+};
+
+
+
template <int dim, int spacedim>
void DataOutBase::write_gnuplot (const typename std::vector<Patch<dim,spacedim> > &patches,
const std::vector<std::string> &data_names,
/* --------------------------- class DataOutInterface ---------------------- */
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::write_dx (std::ostream &out) const
+{
+ DataOutBase::write_dx (get_patches(), get_dataset_names(),
+ dx_flags, out);
+};
+
+
+
template <int dim, int spacedim>
void DataOutInterface<dim,spacedim>::write_ucd (std::ostream &out) const
{
switch (output_format)
{
+ case dx:
+ write_dx (out);
+ break;
+
case ucd:
write_ucd (out);
break;
+template <int dim, int spacedim>
+void DataOutInterface<dim,spacedim>::set_flags (const DXFlags &flags)
+{
+ dx_flags = flags;
+};
+
+
+
template <int dim, int spacedim>
void DataOutInterface<dim,spacedim>::set_flags (const UcdFlags &flags)
{
switch (output_format)
{
+ case dx:
+ return ".dx";
+
case ucd:
return ".inp";
typename DataOutInterface<dim,spacedim>::OutputFormat
DataOutInterface<dim,spacedim>::parse_output_format (const std::string &format_name)
{
+ if (format_name == "dx")
+ return dx;
+
if (format_name == "ucd")
return ucd;
template <int dim, int spacedim>
std::string DataOutInterface<dim,spacedim>::get_output_format_names ()
{
- return "ucd|gnuplot|povray|eps|gmv|vtk";
+ return "dx|ucd|gnuplot|povray|eps|gmv|vtk";
}
prm.declare_entry ("Output format", "gnuplot",
Patterns::Selection (get_output_format_names ()));
+ prm.enter_subsection ("DX output parameters");
+ DXFlags::declare_parameters (prm);
+ prm.leave_subsection ();
+
prm.enter_subsection ("UCD output parameters");
UcdFlags::declare_parameters (prm);
prm.leave_subsection ();
const std::string& output_name = prm.get ("Output format");
default_fmt = parse_output_format (output_name);
+ prm.enter_subsection ("DX output parameters");
+ dx_flags.parse_parameters (prm);
+ prm.leave_subsection ();
+
prm.enter_subsection ("UCD output parameters");
ucd_flags.parse_parameters (prm);
prm.leave_subsection ();
DataOutInterface<dim,spacedim>::memory_consumption () const
{
return (sizeof (default_fmt) +
+ MemoryConsumption::memory_consumption (dx_flags) +
MemoryConsumption::memory_consumption (ucd_flags) +
MemoryConsumption::memory_consumption (gnuplot_flags) +
MemoryConsumption::memory_consumption (povray_flags) +