void initialize_matrices(const std::vector<MatrixPtr>& matrices,
bool both);
+ /**
+ * Initialize quadrature values
+ * to <tt>nv</tt> values in
+ * <tt>np</tt> quadrature points.
+ */
+ void initialize_quadrature(unsigned int np, unsigned int nv);
+
/**
* Reinitialize matrices for
* new cell. Resizes the
*/
unsigned int n_matrices () const;
+ /**
+ * The number of quadrature
+ * points in #quadrature_values.
+ */
+ unsigned int n_quadrature_points() const;
+
+ /**
+ * The number of values in each
+ * quadrature point in
+ * #quadrature_values.
+ */
+ unsigned int n_quadrature_values() const;
+
/**
* Access scalar value at index
* @p i.
*/
const MatrixBlock<FullMatrix<number> >& matrix(unsigned int i, bool external = false) const;
+ /**
+ * Access the <i>i</i>th value
+ * at quadrature point <i>k</i>
+ */
+ number& quadrature_value(unsigned int k, unsigned int i);
+
+ /**
+ * Read the <i>i</i>th value
+ * at quadrature point <i>k</i>
+ */
+ number quadrature_value(unsigned int k, unsigned int i) const;
+
private:
/**
* The local numbers,
* faces.
*/
std::vector<MatrixBlock<FullMatrix<number> > > M2;
+
+ /**
+ * Values in quadrature points.
+ */
+ std::vector<std::vector<number> > quadrature_values;
};
template <int dim, int spacedim> class DoFInfoBox;
* @ingroup MeshWorker
* @author Guido Kanschat, 2009
*/
- template<int dim, int spacedim = dim>
- class DoFInfo : public LocalResults<double>
+ template<int dim, int spacedim = dim, typename number = double>
+ class DoFInfo : public LocalResults<number>
{
public:
/// The current cell
}
}
}
-
+
+ template <typename number>
+ inline void
+ LocalResults<number>::initialize_quadrature(unsigned int np, unsigned int nv)
+ {
+ quadrature_values.resize(np, std::vector<number>(nv));
+ }
+
+
template <typename number>
inline void
LocalResults<number>::reinit(const BlockIndices& bi)
}
+ template <typename number>
+ inline
+ unsigned int
+ LocalResults<number>::n_quadrature_points() const
+ {
+ return quadrature_values.size();
+ }
+
+
+ template <typename number>
+ inline
+ unsigned int
+ LocalResults<number>::n_quadrature_values() const
+ {
+ Assert(quadrature_values.size() != 0, ExcNotInitialized());
+ return quadrature_values[0].size();
+ }
+
+
template <typename number>
inline
number&
return M1[i];
}
+
+ template <typename number>
+ inline
+ number&
+ LocalResults<number>::quadrature_value(unsigned int k, unsigned int i)
+ {
+ AssertIndexRange(k,quadrature_values.size());
+ AssertIndexRange(i,quadrature_values[0].size());
+ return quadrature_values[k][i];
+ }
+
+
template <typename number>
inline
number
return M1[i];
}
+
+ template <typename number>
+ inline
+ number
+ LocalResults<number>::quadrature_value(unsigned int k, unsigned int i) const
+ {
+ AssertIndexRange(k,quadrature_values.size());
+ AssertIndexRange(i,quadrature_values[0].size());
+ return quadrature_values[k][i];
+ }
+
+
//----------------------------------------------------------------------//
- template <int dim, int spacedim>
+ template <int dim, int spacedim, typename number>
template <class DH>
- DoFInfo<dim,spacedim>::DoFInfo(const DH& dof_handler)
+ DoFInfo<dim,spacedim,number>::DoFInfo(const DH& dof_handler)
{
std::vector<unsigned int> aux(1);
aux[0] = dof_handler.get_fe().dofs_per_cell;
}
- template <int dim, int spacedim>
+ template <int dim, int spacedim, typename number>
template <class DHCellIterator>
inline void
- DoFInfo<dim,spacedim>::reinit(const DHCellIterator& c)
+ DoFInfo<dim,spacedim,number>::reinit(const DHCellIterator& c)
{
get_indices(c);
cell = static_cast<typename Triangulation<dim,spacedim>::cell_iterator> (c);
face_number = deal_II_numbers::invalid_unsigned_int;
sub_number = deal_II_numbers::invalid_unsigned_int;
if (block_info)
- LocalResults<double>::reinit(block_info->local());
+ LocalResults<number>::reinit(block_info->local());
else
- LocalResults<double>::reinit(aux_local_indices);
+ LocalResults<number>::reinit(aux_local_indices);
}
- template<int dim, int spacedim>
+ template<int dim, int spacedim, typename number>
template <class DHCellIterator, class DHFaceIterator>
inline void
- DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c,
- const DHFaceIterator& f,
- unsigned int n)
+ DoFInfo<dim,spacedim,number>::reinit(const DHCellIterator& c,
+ const DHFaceIterator& f,
+ unsigned int n)
{
if ((cell.state() != IteratorState::valid)
|| cell != static_cast<typename Triangulation<dim>::cell_iterator> (c))
face_number = n;
sub_number = deal_II_numbers::invalid_unsigned_int;
if (block_info)
- LocalResults<double>::reinit(block_info->local());
+ LocalResults<number>::reinit(block_info->local());
else
- LocalResults<double>::reinit(aux_local_indices);
+ LocalResults<number>::reinit(aux_local_indices);
}
- template<int dim, int spacedim>
+ template<int dim, int spacedim, typename number>
template <class DHCellIterator, class DHFaceIterator>
inline void
- DoFInfo<dim, spacedim>::reinit(const DHCellIterator& c,
- const DHFaceIterator& f,
- unsigned int n,
- const unsigned int s)
+ DoFInfo<dim,spacedim,number>::reinit(const DHCellIterator& c,
+ const DHFaceIterator& f,
+ unsigned int n,
+ const unsigned int s)
{
if (cell.state() != IteratorState::valid
|| cell != static_cast<typename Triangulation<dim>::cell_iterator> (c))
face_number = n;
sub_number = s;
if (block_info)
- LocalResults<double>::reinit(block_info->local());
+ LocalResults<number>::reinit(block_info->local());
else
- LocalResults<double>::reinit(aux_local_indices);
+ LocalResults<number>::reinit(aux_local_indices);
}
- template<int dim, int spacedim>
+ template<int dim, int spacedim, typename number>
inline const BlockIndices&
- DoFInfo<dim, spacedim>::local_indices() const
+ DoFInfo<dim,spacedim,number>::local_indices() const
{
if (block_info)
return block_info->local();
namespace MeshWorker
{
- template <int dim, int spacedim>
- DoFInfo<dim,spacedim>::DoFInfo(const BlockInfo& info)
+ template <int dim, int spacedim, typename number>
+ DoFInfo<dim,spacedim,number>::DoFInfo(const BlockInfo& info)
: block_info(&info, typeid(*this).name())
{}
- template <int dim, int spacedim>
- DoFInfo<dim,spacedim>::DoFInfo()
+ template <int dim, int spacedim, typename number>
+ DoFInfo<dim,spacedim,number>::DoFInfo()
{}
- template <int dim, int spacedim>
+ template <int dim, int spacedim, typename number>
void
- DoFInfo<dim,spacedim>::get_indices(const typename DoFHandler<dim, spacedim>::cell_iterator& c)
+ DoFInfo<dim,spacedim,number>::get_indices(const typename DoFHandler<dim, spacedim>::cell_iterator& c)
{
if (!c->has_children())
{
}
- template <int dim, int spacedim>
+ template <int dim, int spacedim, typename number>
void
- DoFInfo<dim,spacedim>::get_indices(const typename MGDoFHandler<dim, spacedim>::cell_iterator& c)
+ DoFInfo<dim,spacedim,number>::get_indices(const typename MGDoFHandler<dim, spacedim>::cell_iterator& c)
{
indices.resize(c->get_fe().dofs_per_cell);
--- /dev/null
+//---------------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2010 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//---------------------------------------------------------------------------
+
+#ifndef __deal2__mesh_worker_output_h
+#define __deal2__mesh_worker_output_h
+
+#include <numerics/mesh_worker_info.h>
+#include <base/named_data.h>
+#include <base/smartpointer.h>
+#include <lac/block_vector.h>
+#include <multigrid/mg_level_object.h>
+
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace MeshWorker
+{
+ namespace Assembler
+ {
+
+/**
+ * A class that, instead of assembling into a matrix or vector,
+ * outputs the results on a cell to a gnuplot patch. This works only
+ * for elements with support points. The first dim data vectors will
+ * be the coordinates, the following are the data.
+ *
+ * @note In the current implementation, only cell data can be written.
+ */
+ class GnuplotPatch
+ {
+ /**
+ * Constructor.
+ */
+ GnuplotPatch();
+
+ /**
+ * Initialize for writing
+ * <i>n</i> data vectors. The
+ * total number of data
+ * vectors produced is n+dim
+ * and the first dim will be
+ * the space coordinates of
+ * the points.
+ */
+ void initialize(unsigned int dim, unsigned int n_vectors, unsigned int n_points);
+
+ /**
+ * Set the stream #os to
+ * which data is written. If
+ * no stream is selected with
+ * this function, data goes
+ * to #deallog.
+ */
+ void initialize_stream(std::ostream& stream);
+
+ /**
+ * Initialize the local data
+ * in the
+ * DoFInfo
+ * object used later for
+ * assembling.
+ *
+ * The second parameter is
+ * used to distinguish
+ * between the data used on
+ * cells and boundary faces
+ * on the one hand and
+ * interior faces on the
+ * other. Interior faces may
+ * require additional data
+ * being initialized.
+ */
+ template <int dim>
+ void initialize_info(DoFInfo<dim>& info,
+ bool interior_face);
+
+ /**
+ * Assemble the local values
+ * into the global vectors.
+ */
+ template<int dim>
+ void assemble(const DoFInfo<dim>& info);
+
+ /**
+ * Assemble both local values
+ * into the global vectors.
+ */
+ template<int dim>
+ void assemble(const DoFInfo<dim>& info1,
+ const DoFInfo<dim>& info2);
+
+ /**
+ * The value of the ith entry
+ * in #results.
+ */
+ number operator() (unsigned int i) const;
+ private:
+ /**
+ * Write the object T either
+ * to the stream #os, if the
+ * pointer is nonzero, or to
+ * #deallog if no pointer has
+ * been set.
+ */
+ template<typename T>
+ void write(const T& t) const;
+
+ unsigned int patch_dimension;
+ unsigned int n_vectors;
+ unsigned int n_points;
+
+ std::ostream* os;
+ };
+
+//----------------------------------------------------------------------//
+
+ template<typename T>
+ inline void
+ GnuplotPatch::write(const T& d) const
+ {
+ if (os == 0)
+ deallog << d;
+ else
+ (*os) << d;
+ }
+
+
+ inline void
+ GnuplotPatch::GnuplotPatch()
+ :
+ os(0)
+ {}
+
+
+ inline void
+ GnuplotPatch::initialize(unsigned int dim, unsigned int nv, unsigned int np)
+ {
+ patch_dimension = dim;
+ n_vectors = n;
+ n_points = np;
+ }
+
+
+ inline void
+ GnuplotPatch::initialize_stream(std::ostream& stream)
+ {
+ os = &stream;
+ }
+
+
+ template <int dim>
+ inline void
+ GnuplotPatch::initialize_info(DoFInfo<dim>& info,
+ bool interior_face)
+ {
+ AssertDimension(patch_dimension,dim);
+ info.initialize_quadrature(n_points, n_vectors);
+ }
+
+
+ template <int dim>
+ inline void
+ GnuplotPatch::assemble(const DoFInfo<dim>& info)
+ {
+
+ }
+
+
+ template <int dim>
+ inline void
+ GnuplotPatch::assemble(const DoFInfo<dim>& info1, const DoFInfo<dim>& info2)
+ {
+
+ }
+ }
+}
+
+#endif