*/
namespace MeshWorker
{
-/**
- * Worker object for integration of functionals, residuals or matrices.
- *
- * In order to allow for sufficient generality, a few steps have to be
- * undertaken to use this class. The constructor will only fill some
- * default values.
- *
- * First, you should consider if you need values from any vectors in a
- * NamedData object. If so, fill the VectorSelector objects
- * #cell_selector, #boundary_selector and #face_selector with their names
- * and the data type (value, gradient, Hessian) to be extracted.
- *
- * Afterwards, you will need to consider UpdateFlags for FEValues
- * objects. A good start is initialize_update_flags(), which looks at
- * the selectors filled before and adds all the flags needed to get
- * the selection. Additional flags can be set with add_update_flags().
- *
- * Finally, we need to choose quadrature formulas. If you choose to
- * use Gauss formulas only, use initialize_gauss_quadrature() with
- * appropriate values. Otherwise, you can fill the variables
- * #cell_quadrature, #boundary_quadrature and #face_quadrature directly.
- *
- * In order to save time, you can set the variables boundary_fluxes
- * and interior_fluxes of the base class to false, thus telling the
- * Meshworker::loop() not to loop over those faces.
- *
- * All the information in here is used to set up IntegrationInfo
- * objects correctly, typically in an IntegrationInfoBox.
- *
- * @ingroup MeshWorker
- * @author Guido Kanschat, 2009
- */
- template <int dim>
- class IntegrationWorker
- {
- public:
- /**
- * The info type expected by a
- * cell integrator.
- */
- typedef IntegrationInfo<dim, dim> CellInfo;
-
- /**
- * Initialize default values.
- */
- IntegrationWorker();
- /**
- * Initialize the
- * VectorSelector objects
- * #cell_selector,
- * #boundary_selector and
- * #face_selector in order to
- * save computational
- * eeffort. If no selectors
- * are used, then values for
- * all named vectors in
- * DoFInfo::global_data will be
- * computed in all quadrature
- * points.
- *
- * This function will also
- * add UpdateFlags to the
- * flags stored in this class.
- */
- void initialize_update_flags();
-
- /**
- * Add additional values for update.
- */
- void add_update_flags(const UpdateFlags flags, bool cell = true,
- bool boundary = true, bool face = true,
- bool neighbor = true);
-
- /** Assign n-point Gauss
- * quadratures to each of the
- * quadrature rules. Here, a
- * size of zero points means
- * that no loop over these grid
- * entities should be
- * performed.
- */
- void initialize_gauss_quadrature(unsigned int n_cell_points,
- unsigned int n_boundary_points,
- unsigned int n_face_points);
-
- /**
- * Select the vectors from
- * DoFInfo::global_data
- * that should be computed in
- * the quadrature points on cells.
- */
- MeshWorker::VectorSelector cell_selector;
-
- /**
- * Select the vectors from
- * DoFInfo::global_data
- * that should be computed in
- * the quadrature points on
- * boundary faces.
- */
- MeshWorker::VectorSelector boundary_selector;
-
- /**
- * Select the vectors from
- * DoFInfo::global_data
- * that should be computed in
- * the quadrature points on
- * interior faces.
- */
- MeshWorker::VectorSelector face_selector;
-
- /**
- * The set of update flags
- * for boundary cell integration.
- *
- * Defaults to
- * #update_JxW_values.
- */
- UpdateFlags cell_flags;
- /**
- * The set of update flags
- * for boundary face integration.
- *
- * Defaults to
- * #update_JxW_values and
- * #update_normal_vectors.
- */
- UpdateFlags boundary_flags;
-
- /**
- * The set of update flags
- * for interior face integration.
- *
- * Defaults to
- * #update_JxW_values and
- * #update_normal_vectors.
- */
- UpdateFlags face_flags;
-
- /**
- * The set of update flags
- * for interior face integration.
- *
- * Defaults to
- * #update_default, since
- * quadrature weights are
- * taken from the other cell.
- */
- UpdateFlags neighbor_flags;
-
- /**
- * The quadrature rule used
- * on cells.
- */
- Quadrature<dim> cell_quadrature;
-
- /**
- * The quadrature rule used
- * on boundary faces.
- */
- Quadrature<dim-1> boundary_quadrature;
-
- /**
- * The quadrature rule used
- * on interior faces.
- */
- Quadrature<dim-1> face_quadrature;
- };
}
* functions at quadrature points.
* </ol>
*
+ * In order to allow for sufficient generality, a few steps have to be
+ * undertaken to use this class.\18
+ *
+ * First, you should consider if you need values from any vectors in a
+ * NamedData object. If so, fill the VectorSelector objects
+ * #cell_selector, #boundary_selector and #face_selector with their names
+ * and the data type (value, gradient, Hessian) to be extracted.
+ *
+ * Afterwards, you will need to consider UpdateFlags for FEValues
+ * objects. A good start is initialize_update_flags(), which looks at
+ * the selectors filled before and adds all the flags needed to get
+ * the selection. Additional flags can be set with add_update_flags().
+ *
+ * Finally, we need to choose quadrature formulas. If you choose to
+ * use Gauss formulas only, use initialize_gauss_quadrature() with
+ * appropriate values. Otherwise, you can fill the variables
+ * #cell_quadrature, #boundary_quadrature and #face_quadrature directly.
+ *
+ * In order to save time, you can set the variables boundary_fluxes
+ * and interior_fluxes of the base class to false, thus telling the
+ * Meshworker::loop() not to loop over those faces.
+ *
+ * All the information in here is used to set up IntegrationInfo
+ * objects correctly, typically in an IntegrationInfoBox.
+ *
* @ingroup MeshWorker
* @author Guido Kanschat, 2009
*/
/// The type of the info object for cells
typedef IntegrationInfo<dim, spacedim> CellInfo;
- /**
- * @name FEValues setup
- */
- /* @{ */
- template <class WORKER>
- void initialize(const WORKER&,
- const FiniteElement<dim, spacedim>& el,
+ void initialize(const FiniteElement<dim, spacedim>& el,
const Mapping<dim, spacedim>& mapping,
const BlockInfo* block_info = 0);
- template <class WORKER, typename VECTOR>
- void initialize(const WORKER&,
- const FiniteElement<dim, spacedim>& el,
+ template <typename VECTOR>
+ void initialize(const FiniteElement<dim, spacedim>& el,
const Mapping<dim, spacedim>& mapping,
const NamedData<VECTOR*>& data,
const BlockInfo* block_info = 0);
- template <class WORKER, typename VECTOR>
- void initialize(const WORKER&,
- const FiniteElement<dim, spacedim>& el,
+ template <typename VECTOR>
+ void initialize(const FiniteElement<dim, spacedim>& el,
const Mapping<dim, spacedim>& mapping,
const NamedData<MGLevelObject<VECTOR>*>& data,
const BlockInfo* block_info = 0);
+ /**
+ * @name FEValues setup
+ */
+ /* @{ */
+ void initialize_update_flags();
+
+ /**
+ * Add additional values for update.
+ */
+ void add_update_flags(const UpdateFlags flags, bool cell = true,
+ bool boundary = true, bool face = true,
+ bool neighbor = true);
+
+ /** Assign n-point Gauss
+ * quadratures to each of the
+ * quadrature rules. Here, a
+ * size of zero points means
+ * that no loop over these grid
+ * entities should be
+ * performed.
+ */
+ void initialize_gauss_quadrature(unsigned int n_cell_points,
+ unsigned int n_boundary_points,
+ unsigned int n_face_points);
+
+ /**
+ * The set of update flags
+ * for boundary cell integration.
+ *
+ * Defaults to
+ * #update_JxW_values.
+ */
+ UpdateFlags cell_flags;
+ /**
+ * The set of update flags
+ * for boundary face integration.
+ *
+ * Defaults to
+ * #update_JxW_values and
+ * #update_normal_vectors.
+ */
+ UpdateFlags boundary_flags;
+
+ /**
+ * The set of update flags
+ * for interior face integration.
+ *
+ * Defaults to
+ * #update_JxW_values and
+ * #update_normal_vectors.
+ */
+ UpdateFlags face_flags;
+
+ /**
+ * The set of update flags
+ * for interior face integration.
+ *
+ * Defaults to
+ * #update_default, since
+ * quadrature weights are
+ * taken from the other cell.
+ */
+ UpdateFlags neighbor_flags;
+
+ /**
+ * The quadrature rule used
+ * on cells.
+ */
+ Quadrature<dim> cell_quadrature;
+
+ /**
+ * The quadrature rule used
+ * on boundary faces.
+ */
+ Quadrature<dim-1> boundary_quadrature;
+
+ /**
+ * The quadrature rule used
+ * on interior faces.
+ */
+ Quadrature<dim-1> face_quadrature;
/* @} */
/**
*/
/* @{ */
+ /**
+ * Initialize the
+ * VectorSelector objects
+ * #cell_selector,
+ * #boundary_selector and
+ * #face_selector in order to
+ * save computational
+ * eeffort. If no selectors
+ * are used, then values for
+ * all named vectors in
+ * DoFInfo::global_data will be
+ * computed in all quadrature
+ * points.
+ *
+ * This function will also
+ * add UpdateFlags to the
+ * flags stored in this class.
+ */
+ /**
+ * Select the vectors from
+ * DoFInfo::global_data
+ * that should be computed in
+ * the quadrature points on cells.
+ */
+ MeshWorker::VectorSelector cell_selector;
+
+ /**
+ * Select the vectors from
+ * DoFInfo::global_data
+ * that should be computed in
+ * the quadrature points on
+ * boundary faces.
+ */
+ MeshWorker::VectorSelector boundary_selector;
+
+ /**
+ * Select the vectors from
+ * DoFInfo::global_data
+ * that should be computed in
+ * the quadrature points on
+ * interior faces.
+ */
+ MeshWorker::VectorSelector face_selector;
+
boost::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > cell_data;
boost::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > boundary_data;
boost::shared_ptr<MeshWorker::VectorDataBase<dim, spacedim> > face_data;
//----------------------------------------------------------------------//
template <int dim, int sdim>
- template <class WORKER>
- void
- IntegrationInfoBox<dim,sdim>::
- initialize(const WORKER& integrator,
- const FiniteElement<dim,sdim>& el,
- const Mapping<dim,sdim>& mapping,
- const BlockInfo* block_info)
- {
- cell.initialize<FEValues<dim,sdim> >(el, mapping, integrator.cell_quadrature,
- integrator.cell_flags, block_info);
- boundary.initialize<FEFaceValues<dim,sdim> >(el, mapping, integrator.boundary_quadrature,
- integrator.boundary_flags, block_info);
- face.initialize<FEFaceValues<dim,sdim> >(el, mapping, integrator.face_quadrature,
- integrator.face_flags, block_info);
- subface.initialize<FESubfaceValues<dim,sdim> >(el, mapping, integrator.face_quadrature,
- integrator.face_flags, block_info);
- neighbor.initialize<FEFaceValues<dim,sdim> >(el, mapping, integrator.face_quadrature,
- integrator.neighbor_flags, block_info);
- }
-
-
- template <int dim, int sdim>
- template <class WORKER, typename VECTOR>
+ template <typename VECTOR>
void
IntegrationInfoBox<dim,sdim>::initialize(
- const WORKER& integrator,
const FiniteElement<dim,sdim>& el,
const Mapping<dim,sdim>& mapping,
const NamedData<VECTOR*>& data,
const BlockInfo* block_info)
{
- initialize(integrator, el, mapping, block_info);
+ initialize(el, mapping, block_info);
boost::shared_ptr<VectorData<VECTOR, dim, sdim> > p;
- p = boost::shared_ptr<VectorData<VECTOR, dim, sdim> >(new VectorData<VECTOR, dim, sdim> (integrator.cell_selector));
+ p = boost::shared_ptr<VectorData<VECTOR, dim, sdim> >(new VectorData<VECTOR, dim, sdim> (cell_selector));
p->initialize(data);
cell_data = p;
cell.initialize_data(p);
- p = boost::shared_ptr<VectorData<VECTOR, dim, sdim> >(new VectorData<VECTOR, dim, sdim> (integrator.boundary_selector));
+ p = boost::shared_ptr<VectorData<VECTOR, dim, sdim> >(new VectorData<VECTOR, dim, sdim> (boundary_selector));
p->initialize(data);
boundary_data = p;
boundary.initialize_data(p);
- p = boost::shared_ptr<VectorData<VECTOR, dim, sdim> >(new VectorData<VECTOR, dim, sdim> (integrator.face_selector));
+ p = boost::shared_ptr<VectorData<VECTOR, dim, sdim> >(new VectorData<VECTOR, dim, sdim> (face_selector));
p->initialize(data);
face_data = p;
face.initialize_data(p);
template <int dim, int sdim>
- template <class WORKER, typename VECTOR>
+ template <typename VECTOR>
void
IntegrationInfoBox<dim,sdim>::initialize(
- const WORKER& integrator,
const FiniteElement<dim,sdim>& el,
const Mapping<dim,sdim>& mapping,
const NamedData<MGLevelObject<VECTOR>*>& data,
const BlockInfo* block_info)
{
- initialize(integrator, el, mapping, block_info);
+ initialize(el, mapping, block_info);
boost::shared_ptr<MGVectorData<VECTOR, dim, sdim> > p;
- p = boost::shared_ptr<MGVectorData<VECTOR, dim, sdim> >(new MGVectorData<VECTOR, dim, sdim> (integrator.cell_selector));
+ p = boost::shared_ptr<MGVectorData<VECTOR, dim, sdim> >(new MGVectorData<VECTOR, dim, sdim> (cell_selector));
p->initialize(data);
cell_data = p;
cell.initialize_data(p);
- p = boost::shared_ptr<MGVectorData<VECTOR, dim, sdim> >(new MGVectorData<VECTOR, dim, sdim> (integrator.boundary_selector));
+ p = boost::shared_ptr<MGVectorData<VECTOR, dim, sdim> >(new MGVectorData<VECTOR, dim, sdim> (boundary_selector));
p->initialize(data);
boundary_data = p;
boundary.initialize_data(p);
- p = boost::shared_ptr<MGVectorData<VECTOR, dim, sdim> >(new MGVectorData<VECTOR, dim, sdim> (integrator.face_selector));
+ p = boost::shared_ptr<MGVectorData<VECTOR, dim, sdim> >(new MGVectorData<VECTOR, dim, sdim> (face_selector));
p->initialize(data);
face_data = p;
face.initialize_data(p);
subface.initialize_data(p);
neighbor.initialize_data(p);
}
-
-
+
+
template <int dim, int sdim>
template <class DOFINFO>
void
IntegrationInfoBox<dim,sdim>::post_cell(const DoFInfoBox<dim, DOFINFO>&)
{}
-
+
template <int dim, int sdim>
template <class DOFINFO>
void
IntegrationInfoBox<dim,sdim>::post_faces(const DoFInfoBox<dim, DOFINFO>&)
{}
-
-
}
DEAL_II_NAMESPACE_CLOSE
#include <numerics/mesh_worker_info.h>
+#include <base/quadrature_lib.h>
DEAL_II_NAMESPACE_OPEN
0, n_comp, 0, info.indices.size());
}
}
+
+
+//----------------------------------------------------------------------//
+
+
+ template<int dim, int sdim>
+ void
+ IntegrationInfoBox<dim,sdim>::initialize_update_flags ()
+ {
+ cell_flags = update_JxW_values;
+ boundary_flags = UpdateFlags(update_JxW_values | update_normal_vectors);
+ face_flags = boundary_flags;
+ neighbor_flags = update_default;
+
+ if (cell_selector.has_values() != 0) cell_flags |= update_values;
+ if (cell_selector.has_gradients() != 0) cell_flags |= update_gradients;
+ if (cell_selector.has_hessians() != 0) cell_flags |= update_hessians;
+
+ if (boundary_selector.has_values() != 0) boundary_flags |= update_values;
+ if (boundary_selector.has_gradients() != 0) boundary_flags |= update_gradients;
+ if (boundary_selector.has_hessians() != 0) boundary_flags |= update_hessians;
+
+ if (face_selector.has_values() != 0) face_flags |= update_values;
+ if (face_selector.has_gradients() != 0) face_flags |= update_gradients;
+ if (face_selector.has_hessians() != 0) face_flags |= update_hessians;
+
+ if (face_selector.has_values() != 0) neighbor_flags |= update_values;
+ if (face_selector.has_gradients() != 0) neighbor_flags |= update_gradients;
+ if (face_selector.has_hessians() != 0) neighbor_flags |= update_hessians;
+ }
+
+
+ template <int dim, int sdim>
+ void
+ IntegrationInfoBox<dim,sdim>::add_update_flags(
+ const UpdateFlags flags,
+ bool cell,
+ bool boundary,
+ bool face,
+ bool neighbor)
+ {
+ if (cell) cell_flags |= flags;
+ if (boundary) boundary_flags |= flags;
+ if (face) face_flags |= flags;
+ if (neighbor) neighbor_flags |= flags;
+ }
+
+
+ template <int dim, int sdim>
+ void
+ IntegrationInfoBox<dim,sdim>::initialize_gauss_quadrature(
+ unsigned int cp,
+ unsigned int bp,
+ unsigned int fp)
+ {
+ cell_quadrature = QGauss<dim>(cp);
+ boundary_quadrature = QGauss<dim-1>(bp);
+ face_quadrature = QGauss<dim-1>(fp);
+
+ }
+
+
+ template <int dim, int sdim>
+ void
+ IntegrationInfoBox<dim,sdim>::
+ initialize(const FiniteElement<dim,sdim>& el,
+ const Mapping<dim,sdim>& mapping,
+ const BlockInfo* block_info)
+ {
+ cell.initialize<FEValues<dim,sdim> >(el, mapping, cell_quadrature,
+ cell_flags, block_info);
+ boundary.initialize<FEFaceValues<dim,sdim> >(el, mapping, boundary_quadrature,
+ boundary_flags, block_info);
+ face.initialize<FEFaceValues<dim,sdim> >(el, mapping, face_quadrature,
+ face_flags, block_info);
+ subface.initialize<FESubfaceValues<dim,sdim> >(el, mapping, face_quadrature,
+ face_flags, block_info);
+ neighbor.initialize<FEFaceValues<dim,sdim> >(el, mapping, face_quadrature,
+ neighbor_flags, block_info);
+ }
}
DEAL_II_NAMESPACE_OPEN
-namespace MeshWorker
-{
-
- template <int dim>
- IntegrationWorker<dim>::IntegrationWorker ()
- {
- cell_flags = update_JxW_values;
- boundary_flags = UpdateFlags(update_JxW_values | update_normal_vectors);
- face_flags = boundary_flags;
- neighbor_flags = update_default;
- }
-
-
- template<int dim>
- void
- IntegrationWorker<dim>::initialize_update_flags ()
- {
- if (cell_selector.has_values() != 0) cell_flags |= update_values;
- if (cell_selector.has_gradients() != 0) cell_flags |= update_gradients;
- if (cell_selector.has_hessians() != 0) cell_flags |= update_hessians;
-
- if (boundary_selector.has_values() != 0) boundary_flags |= update_values;
- if (boundary_selector.has_gradients() != 0) boundary_flags |= update_gradients;
- if (boundary_selector.has_hessians() != 0) boundary_flags |= update_hessians;
-
- if (face_selector.has_values() != 0) face_flags |= update_values;
- if (face_selector.has_gradients() != 0) face_flags |= update_gradients;
- if (face_selector.has_hessians() != 0) face_flags |= update_hessians;
-
- if (face_selector.has_values() != 0) neighbor_flags |= update_values;
- if (face_selector.has_gradients() != 0) neighbor_flags |= update_gradients;
- if (face_selector.has_hessians() != 0) neighbor_flags |= update_hessians;
- }
-
-
- template<int dim>
- void
- IntegrationWorker<dim>::add_update_flags(
- const UpdateFlags flags, bool cell, bool boundary, bool face, bool neighbor)
- {
- if (cell) cell_flags |= flags;
- if (boundary) boundary_flags |= flags;
- if (face) face_flags |= flags;
- if (neighbor) neighbor_flags |= flags;
- }
-
-
- template<int dim>
- void
- IntegrationWorker<dim>::initialize_gauss_quadrature(
- unsigned int cp,
- unsigned int bp,
- unsigned int fp)
- {
- cell_quadrature = QGauss<dim>(cp);
- boundary_quadrature = QGauss<dim-1>(bp);
- face_quadrature = QGauss<dim-1>(fp);
-
- }
-}
-
-#if deal_II_dimension > 1
-namespace MeshWorker
-{
- template class IntegrationWorker<deal_II_dimension>;
-}
-#endif
-
DEAL_II_NAMESPACE_CLOSE
namespace MeshWorker
{
template class IntegrationInfo<deal_II_dimension, deal_II_dimension>;
+ template class IntegrationInfoBox<deal_II_dimension, deal_II_dimension>;
template class LocalResults<float>;
template class DoFInfo<deal_II_dimension,deal_II_dimension,float>;
// order to make our life easier
// below.
typedef MeshWorker::DoFInfo<dim> DoFInfo;
- typedef typename MeshWorker::IntegrationWorker<dim>::CellInfo CellInfo;
+ typedef MeshWorker::IntegrationInfo<dim> CellInfo;
// The following three functions
// are then the ones that get called
// object distributes these into
// the global sparse matrix and the
// right hand side vector.
- MeshWorker::IntegrationWorker<dim> integration_worker;
-
+ MeshWorker::IntegrationInfoBox<dim> info_box;
+
// First, we initialize the
// quadrature formulae and the
// update flags in the worker base
// independently, we have to hand
// over this value three times.
const unsigned int n_gauss_points = dof_handler.get_fe().degree+1;
- integration_worker.initialize_gauss_quadrature(n_gauss_points,
+ info_box.initialize_gauss_quadrature(n_gauss_points,
n_gauss_points,
n_gauss_points);
-
+
// These are the types of values we
// need for integrating our
// system. They are added to the
// flags used on cells, boundary
// and interior faces, as well as
// interior neighbor faces, which is
- // forced by the four @p true values.
+ // forced by the four @p true
+ // values.
+ info_box.initialize_update_flags();
UpdateFlags update_flags = update_quadrature_points |
update_values |
update_gradients;
- integration_worker.add_update_flags(update_flags, true, true, true, true);
+ info_box.add_update_flags(update_flags, true, true, true, true);
+
+ // After preparing all data in
+ // <tt>info_box</tt>, we initialize
+ // the FEValus objects in there.
+ info_box.initialize(fe, mapping);
+
+ // The object created so far helps
+ // us do the local integration on
+ // each cell and face. Now, we need
+ // an object which receives the
+ // integrated (local) data and
+ // forwards them to the assembler.
+ MeshWorker::DoFInfo<dim> dof_info(dof_handler);
- // Finally, we have to tell the
- // assembler base class where to
- // put the local data. These will
- // be our system matrix and the
- // right hand side.
+ // Now, we have to create the
+ // assembler object and tell it,
+ // where to put the local
+ // data. These will be our system
+ // matrix and the right hand side.
MeshWorker::Assembler::SystemSimple<SparseMatrix<double>, Vector<double> >
assembler;
assembler.initialize(system_matrix, right_hand_side);
-
- // We are now ready to get to the
- // integration loop. @p info_box is
- // an object that generates the
- // extended iterators for cells and
- // faces of type
- // MeshWorker::IntegrationInfo. Since
- // we need five different of them,
- // this is a handy shortcut. It
- // receives all the stuff we
- // created so far.
- MeshWorker::IntegrationInfoBox<dim> info_box;
- MeshWorker::DoFInfo<dim> dof_info(dof_handler);
- info_box.initialize(integration_worker, fe, mapping);
-
+
// Finally, the integration loop
// over all active cells
// (determined by the first
// interior faces, respectively. All
// the information needed for the
// local integration is provided by
- // MeshWorker::IntegrationWorker<dim>::CellInfo. Note
+ // MeshWorker::IntegrationInfo<dim>. Note
// that this public interface cannot
// be changed, because it is expected
// by MeshWorker::integration_loop().
{
public:
static void cell(MeshWorker::DoFInfo<dim>& dinfo,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info);
+ typename MeshWorker::IntegrationInfo<dim>& info);
static void bdry(MeshWorker::DoFInfo<dim>& dinfo,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info);
+ typename MeshWorker::IntegrationInfo<dim>& info);
static void face(MeshWorker::DoFInfo<dim>& dinfo1,
MeshWorker::DoFInfo<dim>& dinfo2,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info1,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info2);
+ typename MeshWorker::IntegrationInfo<dim>& info1,
+ typename MeshWorker::IntegrationInfo<dim>& info2);
};
template <int dim>
void MatrixIntegrator<dim>::cell(
MeshWorker::DoFInfo<dim>& dinfo,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info)
+ typename MeshWorker::IntegrationInfo<dim>& info)
{
const FEValuesBase<dim>& fe = info.fe_values();
FullMatrix<double>& local_matrix = dinfo.matrix(0,false).matrix;
template <int dim>
void MatrixIntegrator<dim>::bdry(
MeshWorker::DoFInfo<dim>& dinfo,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info)
+ typename MeshWorker::IntegrationInfo<dim>& info)
{
const FEValuesBase<dim>& fe = info.fe_values();
FullMatrix<double>& local_matrix = dinfo.matrix(0,false).matrix;
void MatrixIntegrator<dim>::face(
MeshWorker::DoFInfo<dim>& dinfo1,
MeshWorker::DoFInfo<dim>& dinfo2,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info1,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info2)
+ typename MeshWorker::IntegrationInfo<dim>& info1,
+ typename MeshWorker::IntegrationInfo<dim>& info2)
{
const FEValuesBase<dim>& fe1 = info1.fe_values();
const FEValuesBase<dim>& fe2 = info2.fe_values();
class RHSIntegrator : public Subscriptor
{
public:
- static void cell(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationWorker<dim>::CellInfo& info);
- static void bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationWorker<dim>::CellInfo& info);
+ static void cell(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationInfo<dim>& info);
+ static void bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationInfo<dim>& info);
static void face(MeshWorker::DoFInfo<dim>& dinfo1,
MeshWorker::DoFInfo<dim>& dinfo2,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info1,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info2);
+ typename MeshWorker::IntegrationInfo<dim>& info1,
+ typename MeshWorker::IntegrationInfo<dim>& info2);
};
template <int dim>
-void RHSIntegrator<dim>::cell(MeshWorker::DoFInfo<dim>&, typename MeshWorker::IntegrationWorker<dim>::CellInfo&)
+void RHSIntegrator<dim>::cell(MeshWorker::DoFInfo<dim>&, typename MeshWorker::IntegrationInfo<dim>&)
{}
template <int dim>
-void RHSIntegrator<dim>::bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationWorker<dim>::CellInfo& info)
+void RHSIntegrator<dim>::bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationInfo<dim>& info)
{
const FEValuesBase<dim>& fe = info.fe_values();
Vector<double>& local_vector = dinfo.vector(0).block(0);
template <int dim>
void RHSIntegrator<dim>::face(MeshWorker::DoFInfo<dim>&,
MeshWorker::DoFInfo<dim>&,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo&,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo&)
+ typename MeshWorker::IntegrationInfo<dim>&,
+ typename MeshWorker::IntegrationInfo<dim>&)
{}
class Estimator : public Subscriptor
{
public:
- static void cell(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationWorker<dim>::CellInfo& info);
- static void bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationWorker<dim>::CellInfo& info);
+ static void cell(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationInfo<dim>& info);
+ static void bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationInfo<dim>& info);
static void face(MeshWorker::DoFInfo<dim>& dinfo1,
MeshWorker::DoFInfo<dim>& dinfo2,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info1,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info2);
+ typename MeshWorker::IntegrationInfo<dim>& info1,
+ typename MeshWorker::IntegrationInfo<dim>& info2);
};
template <int dim>
-void Estimator<dim>::cell(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationWorker<dim>::CellInfo& info)
+void Estimator<dim>::cell(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationInfo<dim>& info)
{
const FEValuesBase<dim>& fe = info.fe_values();
template <int dim>
-void Estimator<dim>::bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationWorker<dim>::CellInfo& info)
+void Estimator<dim>::bdry(MeshWorker::DoFInfo<dim>& dinfo, typename MeshWorker::IntegrationInfo<dim>& info)
{
const FEValuesBase<dim>& fe = info.fe_values();
template <int dim>
void Estimator<dim>::face(MeshWorker::DoFInfo<dim>& dinfo1,
MeshWorker::DoFInfo<dim>& dinfo2,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info1,
- typename MeshWorker::IntegrationWorker<dim>::CellInfo& info2)
+ typename MeshWorker::IntegrationInfo<dim>& info1,
+ typename MeshWorker::IntegrationInfo<dim>& info2)
{
const FEValuesBase<dim>& fe = info1.fe_values();
const std::vector<double>& uh1 = info1.values[0][0];
class Step39
{
public:
- typedef typename MeshWorker::IntegrationWorker<dim>::CellInfo CellInfo;
+ typedef typename MeshWorker::IntegrationInfo<dim> CellInfo;
Step39(const FiniteElement<dim>& fe);
void
Step39<dim>::assemble_matrix()
{
- MeshWorker::IntegrationWorker<dim> integration_worker;
- MeshWorker::Assembler::MatrixSimple<SparseMatrix<double> > assembler;
-
+ MeshWorker::IntegrationInfoBox<dim> info_box;
const unsigned int n_gauss_points = dof_handler.get_fe().tensor_degree()+1;
- integration_worker.initialize_gauss_quadrature(n_gauss_points, n_gauss_points, n_gauss_points);
+ info_box.initialize_gauss_quadrature(n_gauss_points, n_gauss_points, n_gauss_points);
+ info_box.initialize_update_flags();
UpdateFlags update_flags = update_values | update_gradients;
- integration_worker.add_update_flags(update_flags, true, true, true, true);
-
- assembler.initialize(matrix);
- MeshWorker::IntegrationInfoBox<dim> info_box;
+ info_box.add_update_flags(update_flags, true, true, true, true);
+ info_box.initialize(fe, mapping);
+
MeshWorker::DoFInfo<dim> dof_info(dof_handler);
- info_box.initialize(integration_worker, fe, mapping);
+
+ MeshWorker::Assembler::MatrixSimple<SparseMatrix<double> > assembler;
+ assembler.initialize(matrix);
+
MeshWorker::integration_loop<dim, dim>(
dof_handler.begin_active(), dof_handler.end(),
dof_info, info_box,
void
Step39<dim>::assemble_mg_matrix()
{
- MeshWorker::IntegrationWorker<dim> integration_worker;
- MeshWorker::Assembler::MGMatrixSimple<SparseMatrix<double> > assembler;
-
+ MeshWorker::IntegrationInfoBox<dim> info_box;
const unsigned int n_gauss_points = mg_dof_handler.get_fe().tensor_degree()+1;
- integration_worker.initialize_gauss_quadrature(n_gauss_points, n_gauss_points, n_gauss_points);
+ info_box.initialize_gauss_quadrature(n_gauss_points, n_gauss_points, n_gauss_points);
+ info_box.initialize_update_flags();
UpdateFlags update_flags = update_values | update_gradients;
- integration_worker.add_update_flags(update_flags, true, true, true, true);
+ info_box.add_update_flags(update_flags, true, true, true, true);
+ info_box.initialize(fe, mapping);
+ MeshWorker::DoFInfo<dim> dof_info(mg_dof_handler);
+
+ MeshWorker::Assembler::MGMatrixSimple<SparseMatrix<double> > assembler;
assembler.initialize(mg_matrix);
assembler.initialize_fluxes(mg_matrix_dg_up, mg_matrix_dg_down);
- MeshWorker::IntegrationInfoBox<dim> info_box;
- MeshWorker::DoFInfo<dim> dof_info(mg_dof_handler);
- info_box.initialize(integration_worker, fe, mapping);
+
MeshWorker::integration_loop<dim, dim> (
mg_dof_handler.begin(), mg_dof_handler.end(),
dof_info, info_box,
void
Step39<dim>::assemble_right_hand_side()
{
- MeshWorker::IntegrationWorker<dim> integration_worker;
- MeshWorker::Assembler::ResidualSimple<Vector<double> > assembler;
-
+ MeshWorker::IntegrationInfoBox<dim> info_box;
const unsigned int n_gauss_points = dof_handler.get_fe().tensor_degree()+1;
- integration_worker.initialize_gauss_quadrature(n_gauss_points, n_gauss_points+1, n_gauss_points);
+ info_box.initialize_gauss_quadrature(n_gauss_points, n_gauss_points+1, n_gauss_points);
+ info_box.initialize_update_flags();
UpdateFlags update_flags = update_quadrature_points | update_values | update_gradients;
- integration_worker.add_update_flags(update_flags, true, true, true, true);
-
+ info_box.add_update_flags(update_flags, true, true, true, true);
+ info_box.initialize(fe, mapping);
+
+ MeshWorker::DoFInfo<dim> dof_info(dof_handler);
+
+ MeshWorker::Assembler::ResidualSimple<Vector<double> > assembler;
NamedData<Vector<double>* > data;
Vector<double>* rhs = &right_hand_side;
data.add(rhs, "RHS");
assembler.initialize(data);
- MeshWorker::IntegrationInfoBox<dim> info_box;
- MeshWorker::DoFInfo<dim> dof_info(dof_handler);
- info_box.initialize(integration_worker, fe, mapping);
MeshWorker::integration_loop<dim, dim>(
dof_handler.begin_active(), dof_handler.end(),
cell != triangulation.end();++cell,++i)
cell->set_user_index(i);
- MeshWorker::IntegrationWorker<dim> integration_worker;
- MeshWorker::Assembler::CellsAndFaces<double> assembler;
-
+ MeshWorker::IntegrationInfoBox<dim> info_box;
const unsigned int n_gauss_points = dof_handler.get_fe().tensor_degree()+1;
- integration_worker.initialize_gauss_quadrature(n_gauss_points, n_gauss_points+1, n_gauss_points);
+ info_box.initialize_gauss_quadrature(n_gauss_points, n_gauss_points+1, n_gauss_points);
NamedData<Vector<double>* > solution_data;
solution_data.add(&solution, "solution");
cs.add("solution", true, true, true);
fs.add("solution", true, true, false);
- integration_worker.cell_selector = cs;
- integration_worker.boundary_selector = fs;
- integration_worker.face_selector = fs;
- integration_worker.initialize_update_flags();
- integration_worker.add_update_flags(update_quadrature_points, false, true, false, false);
-
+ info_box.cell_selector = cs;
+ info_box.boundary_selector = fs;
+ info_box.face_selector = fs;
+ info_box.initialize_update_flags();
+ info_box.add_update_flags(update_quadrature_points, false, true, false, false);
+ info_box.initialize(fe, mapping, solution_data);
+
+ MeshWorker::DoFInfo<dim> dof_info(dof_handler);
+
+ MeshWorker::Assembler::CellsAndFaces<double> assembler;
NamedData<BlockVector<double>* > out_data;
BlockVector<double>* est = &estimates;
out_data.add(est, "cells");
assembler.initialize(out_data, false);
- MeshWorker::IntegrationInfoBox<dim> info_box;
- MeshWorker::DoFInfo<dim> dof_info(dof_handler);
- info_box.initialize(integration_worker, fe, mapping, solution_data);
+
MeshWorker::integration_loop<dim, dim> (
dof_handler.begin_active(), dof_handler.end(),
dof_info, info_box,