From 0e73b9b3edc9c97fd78c1d57743229bf079607d5 Mon Sep 17 00:00:00 2001 From: kanschat Date: Tue, 16 Mar 2010 15:48:35 +0000 Subject: [PATCH] prepare for output through MeshWorker git-svn-id: https://svn.dealii.org/trunk@20831 0785d39b-7218-0410-832d-ea1e28bc413d --- .../include/numerics/mesh_worker_info.h | 136 ++++++++++--- .../numerics/mesh_worker_info.templates.h | 16 +- .../include/numerics/mesh_worker_output.h | 187 ++++++++++++++++++ .../source/numerics/mesh_worker_info.cc | 5 +- 4 files changed, 311 insertions(+), 33 deletions(-) create mode 100644 deal.II/deal.II/include/numerics/mesh_worker_output.h diff --git a/deal.II/deal.II/include/numerics/mesh_worker_info.h b/deal.II/deal.II/include/numerics/mesh_worker_info.h index f20fe04124..b591f74641 100644 --- a/deal.II/deal.II/include/numerics/mesh_worker_info.h +++ b/deal.II/deal.II/include/numerics/mesh_worker_info.h @@ -88,6 +88,13 @@ namespace MeshWorker void initialize_matrices(const std::vector& matrices, bool both); + /** + * Initialize quadrature values + * to nv values in + * np quadrature points. + */ + void initialize_quadrature(unsigned int np, unsigned int nv); + /** * Reinitialize matrices for * new cell. Resizes the @@ -111,6 +118,19 @@ namespace MeshWorker */ 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. @@ -155,6 +175,18 @@ namespace MeshWorker */ const MatrixBlock >& matrix(unsigned int i, bool external = false) const; + /** + * Access the ith value + * at quadrature point k + */ + number& quadrature_value(unsigned int k, unsigned int i); + + /** + * Read the ith value + * at quadrature point k + */ + number quadrature_value(unsigned int k, unsigned int i) const; + private: /** * The local numbers, @@ -191,6 +223,11 @@ namespace MeshWorker * faces. */ std::vector > > M2; + + /** + * Values in quadrature points. + */ + std::vector > quadrature_values; }; template class DoFInfoBox; @@ -229,8 +266,8 @@ namespace MeshWorker * @ingroup MeshWorker * @author Guido Kanschat, 2009 */ - template - class DoFInfo : public LocalResults + template + class DoFInfo : public LocalResults { public: /// The current cell @@ -788,8 +825,16 @@ namespace MeshWorker } } } - + + template + inline void + LocalResults::initialize_quadrature(unsigned int np, unsigned int nv) + { + quadrature_values.resize(np, std::vector(nv)); + } + + template inline void LocalResults::reinit(const BlockIndices& bi) @@ -834,6 +879,25 @@ namespace MeshWorker } + template + inline + unsigned int + LocalResults::n_quadrature_points() const + { + return quadrature_values.size(); + } + + + template + inline + unsigned int + LocalResults::n_quadrature_values() const + { + Assert(quadrature_values.size() != 0, ExcNotInitialized()); + return quadrature_values[0].size(); + } + + template inline number& @@ -868,6 +932,18 @@ namespace MeshWorker return M1[i]; } + + template + inline + number& + LocalResults::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 inline number @@ -902,11 +978,23 @@ namespace MeshWorker return M1[i]; } + + template + inline + number + LocalResults::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 + template template - DoFInfo::DoFInfo(const DH& dof_handler) + DoFInfo::DoFInfo(const DH& dof_handler) { std::vector aux(1); aux[0] = dof_handler.get_fe().dofs_per_cell; @@ -914,28 +1002,28 @@ namespace MeshWorker } - template + template template inline void - DoFInfo::reinit(const DHCellIterator& c) + DoFInfo::reinit(const DHCellIterator& c) { get_indices(c); cell = static_cast::cell_iterator> (c); face_number = deal_II_numbers::invalid_unsigned_int; sub_number = deal_II_numbers::invalid_unsigned_int; if (block_info) - LocalResults::reinit(block_info->local()); + LocalResults::reinit(block_info->local()); else - LocalResults::reinit(aux_local_indices); + LocalResults::reinit(aux_local_indices); } - template + template template inline void - DoFInfo::reinit(const DHCellIterator& c, - const DHFaceIterator& f, - unsigned int n) + DoFInfo::reinit(const DHCellIterator& c, + const DHFaceIterator& f, + unsigned int n) { if ((cell.state() != IteratorState::valid) || cell != static_cast::cell_iterator> (c)) @@ -945,19 +1033,19 @@ namespace MeshWorker face_number = n; sub_number = deal_II_numbers::invalid_unsigned_int; if (block_info) - LocalResults::reinit(block_info->local()); + LocalResults::reinit(block_info->local()); else - LocalResults::reinit(aux_local_indices); + LocalResults::reinit(aux_local_indices); } - template + template template inline void - DoFInfo::reinit(const DHCellIterator& c, - const DHFaceIterator& f, - unsigned int n, - const unsigned int s) + DoFInfo::reinit(const DHCellIterator& c, + const DHFaceIterator& f, + unsigned int n, + const unsigned int s) { if (cell.state() != IteratorState::valid || cell != static_cast::cell_iterator> (c)) @@ -967,15 +1055,15 @@ namespace MeshWorker face_number = n; sub_number = s; if (block_info) - LocalResults::reinit(block_info->local()); + LocalResults::reinit(block_info->local()); else - LocalResults::reinit(aux_local_indices); + LocalResults::reinit(aux_local_indices); } - template + template inline const BlockIndices& - DoFInfo::local_indices() const + DoFInfo::local_indices() const { if (block_info) return block_info->local(); diff --git a/deal.II/deal.II/include/numerics/mesh_worker_info.templates.h b/deal.II/deal.II/include/numerics/mesh_worker_info.templates.h index 28c1afe4f9..26eccbb9ed 100644 --- a/deal.II/deal.II/include/numerics/mesh_worker_info.templates.h +++ b/deal.II/deal.II/include/numerics/mesh_worker_info.templates.h @@ -19,20 +19,20 @@ DEAL_II_NAMESPACE_OPEN namespace MeshWorker { - template - DoFInfo::DoFInfo(const BlockInfo& info) + template + DoFInfo::DoFInfo(const BlockInfo& info) : block_info(&info, typeid(*this).name()) {} - template - DoFInfo::DoFInfo() + template + DoFInfo::DoFInfo() {} - template + template void - DoFInfo::get_indices(const typename DoFHandler::cell_iterator& c) + DoFInfo::get_indices(const typename DoFHandler::cell_iterator& c) { if (!c->has_children()) { @@ -53,9 +53,9 @@ namespace MeshWorker } - template + template void - DoFInfo::get_indices(const typename MGDoFHandler::cell_iterator& c) + DoFInfo::get_indices(const typename MGDoFHandler::cell_iterator& c) { indices.resize(c->get_fe().dofs_per_cell); diff --git a/deal.II/deal.II/include/numerics/mesh_worker_output.h b/deal.II/deal.II/include/numerics/mesh_worker_output.h new file mode 100644 index 0000000000..68bfcd2fec --- /dev/null +++ b/deal.II/deal.II/include/numerics/mesh_worker_output.h @@ -0,0 +1,187 @@ +//--------------------------------------------------------------------------- +// $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 +#include +#include +#include +#include + + +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 + * n 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 + void initialize_info(DoFInfo& info, + bool interior_face); + + /** + * Assemble the local values + * into the global vectors. + */ + template + void assemble(const DoFInfo& info); + + /** + * Assemble both local values + * into the global vectors. + */ + template + void assemble(const DoFInfo& info1, + const DoFInfo& 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 + void write(const T& t) const; + + unsigned int patch_dimension; + unsigned int n_vectors; + unsigned int n_points; + + std::ostream* os; + }; + +//----------------------------------------------------------------------// + + template + 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 + inline void + GnuplotPatch::initialize_info(DoFInfo& info, + bool interior_face) + { + AssertDimension(patch_dimension,dim); + info.initialize_quadrature(n_points, n_vectors); + } + + + template + inline void + GnuplotPatch::assemble(const DoFInfo& info) + { + + } + + + template + inline void + GnuplotPatch::assemble(const DoFInfo& info1, const DoFInfo& info2) + { + + } + } +} + +#endif diff --git a/deal.II/deal.II/source/numerics/mesh_worker_info.cc b/deal.II/deal.II/source/numerics/mesh_worker_info.cc index 663b29e61b..c53b15956a 100644 --- a/deal.II/deal.II/source/numerics/mesh_worker_info.cc +++ b/deal.II/deal.II/source/numerics/mesh_worker_info.cc @@ -26,8 +26,11 @@ DEAL_II_NAMESPACE_OPEN namespace MeshWorker { + template class LocalResults; + template class DoFInfo; + template class LocalResults; - template class DoFInfo; + template class DoFInfo; template class IntegrationInfo; template void IntegrationInfo -- 2.39.5