--- /dev/null
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2000 - 2013 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// test meshworker LoopControl
+// variation of mesh_worker_01 with more cpus and cells
+
+#include "../tests.h"
+#include <deal.II/meshworker/loop.h>
+#include <deal.II/meshworker/assembler.h>
+
+#include <deal.II/fe/fe_dgp.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/filtered_iterator.h>
+
+#include <fstream>
+#include <iomanip>
+
+using namespace dealii;
+
+
+template <int dim>
+class myIntegrator: public dealii::MeshWorker::LocalIntegrator<dim>
+{
+public:
+ typedef MeshWorker::IntegrationInfo<dim> CellInfo;
+
+ void cell(MeshWorker::DoFInfo<dim> &dinfo, CellInfo &info) const;
+ void boundary(MeshWorker::DoFInfo<dim> &dinfo, CellInfo &info) const;
+ void face(MeshWorker::DoFInfo<dim> &dinfo1, MeshWorker::DoFInfo<dim> &dinfo2,
+ CellInfo &info1, CellInfo &info2) const;
+
+
+};
+
+template <int dim>
+void
+myIntegrator<dim>::cell(MeshWorker::DoFInfo<dim> &info, CellInfo &) const
+{
+ deallog << "C " << info.cell->id() << std::endl;
+}
+
+
+template <int dim>
+void
+myIntegrator<dim>::boundary(MeshWorker::DoFInfo<dim> &info, CellInfo &) const
+{
+ //deallog << "B cell = " << info.cell->id() << " face = " << info.face_number << std::endl;
+}
+
+
+template <int dim>
+void
+myIntegrator<dim>::face(MeshWorker::DoFInfo<dim> &info1, MeshWorker::DoFInfo<dim> &info2,
+ CellInfo &, CellInfo &) const
+{
+ deallog << "F cell1 = " << info1.cell->id()
+ << " face = " << info1.face_number
+ << " cell2 = " << info2.cell->id()
+ << " face2 = " << info2.face_number
+ << std::endl;
+}
+
+
+class DoNothingAssembler
+{
+ public:
+ template <class DOFINFO>
+ void initialize_info(DOFINFO &info, bool face) const {}
+ template<class DOFINFO>
+ void assemble(const DOFINFO &info){}
+ template<class DOFINFO>
+ void assemble(const DOFINFO &info1,
+ const DOFINFO &info2) {}
+
+
+ };
+
+template <int dim>
+void
+test_simple(DoFHandler<dim> &dofs, MeshWorker::LoopControl &lctrl)
+{
+ myIntegrator<dim> local;
+ DoNothingAssembler assembler;
+ MeshWorker::IntegrationInfoBox<dim> info_box;
+
+ MeshWorker::DoFInfo<dim> dof_info(dofs.block_info());
+
+// integration_loop(ITERATOR begin,
+// typename identity<ITERATOR>::type end,
+// DOFINFO &dinfo,
+// INFOBOX &info,
+// const std_cxx1x::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker,
+// const std_cxx1x::function<void (DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker,
+// const std_cxx1x::function<void (DOFINFO &, DOFINFO &,
+// typename INFOBOX::CellInfo &,
+// typename INFOBOX::CellInfo &)> &face_worker,
+// ASSEMBLER &assembler,
+// const LoopControl &lctrl)
+//
+
+
+ MeshWorker::integration_loop<dim, dim, typename DoFHandler<dim>::active_cell_iterator, DoNothingAssembler>
+ (dofs.begin_active(), dofs.end(),
+ dof_info, info_box,
+ local,
+ assembler,
+ lctrl);
+
+// MeshWorker::loop<dim, dim, MeshWorker::DoFInfo<dim>, MeshWorker::IntegrationInfoBox<dim> >
+// (dofs.begin_active(), dofs.end(),
+// dof_info, info_box,
+// std_cxx1x::bind (&Integrator<dim>::cell, local, std_cxx1x::_1, std_cxx1x::_2),
+// std_cxx1x::bind (&Integrator<dim>::bdry, local, std_cxx1x::_1, std_cxx1x::_2),
+// std_cxx1x::bind (&Integrator<dim>::face, local, std_cxx1x::_1, std_cxx1x::_2, std_cxx1x::_3, std_cxx1x::_4),
+// local,
+// lctrl);
+}
+
+std::string id_to_string(const CellId &id)
+{
+ std::ostringstream ss;
+ ss << id;
+ return ss.str();
+}
+
+template<int dim>
+void test_loop(DoFHandler<dim> &dofs, MeshWorker::LoopControl &lctrl)
+{
+ deallog << "* own_cells=" << lctrl.own_cells
+ << " ghost_cells=" << lctrl.ghost_cells
+ << " own_faces=" << lctrl.own_faces
+ << " faces_to_ghost=" << lctrl.faces_to_ghost
+ << std::endl;
+ test_simple(dofs, lctrl);
+}
+
+template<int dim>
+void
+test()
+{
+ parallel::distributed::Triangulation<dim> tr(MPI_COMM_WORLD,
+ Triangulation<dim>::none/*,
+ parallel::distributed::Triangulation<dim>::construct_multigrid_hierarchy*/);
+ GridGenerator::hyper_cube(tr);
+ tr.refine_global(2);
+
+ FE_DGP<dim> fe(0);
+
+ DoFHandler<dim> dofs(tr);
+ dofs.distribute_dofs(fe);
+
+ dofs.initialize_local_block_info();
+ deallog << "DoFHandler ndofs=" << dofs.n_dofs() << std::endl;
+
+ MeshWorker::LoopControl lctrl;
+
+ deallog << "*** 1. CELLS ***" << std::endl;
+ /*
+ lctrl.own_faces = MeshWorker::LoopControl::never;
+ lctrl.faces_to_ghost = MeshWorker::LoopControl::never;
+
+ lctrl.own_cells = false; lctrl.ghost_cells = false;
+ test_loop(dofs, lctrl);
+
+ lctrl.own_cells = true; lctrl.ghost_cells = false;
+ test_loop(dofs, lctrl);
+
+ lctrl.own_cells = false; lctrl.ghost_cells = true;
+ test_loop(dofs, lctrl);
+
+ lctrl.own_cells = true; lctrl.ghost_cells = true;
+ test_loop(dofs, lctrl);
+ */
+ deallog << "*** 2. FACES ***" << std::endl;
+
+ lctrl.own_cells = false; lctrl.ghost_cells = false;
+
+ lctrl.own_faces = MeshWorker::LoopControl::one;
+ lctrl.faces_to_ghost = MeshWorker::LoopControl::never;
+ test_loop(dofs, lctrl);
+
+ lctrl.own_faces = MeshWorker::LoopControl::both;
+ lctrl.faces_to_ghost = MeshWorker::LoopControl::never;
+ test_loop(dofs, lctrl);
+
+ lctrl.own_faces = MeshWorker::LoopControl::never;
+ lctrl.faces_to_ghost = MeshWorker::LoopControl::one;
+ test_loop(dofs, lctrl);
+
+ lctrl.own_faces = MeshWorker::LoopControl::never;
+ lctrl.faces_to_ghost = MeshWorker::LoopControl::both;
+ test_loop(dofs, lctrl);
+
+//
+//
+// for (int gc=0;gc<2;gc++)
+// for (int oc=0;oc<2;oc++)
+// for (int of=0;of<3;of++)
+// {
+// lctrl.own_cells = !!oc;
+// lctrl.ghost_cells = !!gc;
+//
+// lctrl.own_faces = (MeshWorker::LoopControl::FaceOption)of;
+// test_loop(dofs, lctrl);
+// }
+}
+
+
+int main (int argc, char **argv)
+{
+ Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv);
+ MPILogInitAll log;
+
+ test<2>();
+}
--- /dev/null
+
+DEAL:0::DoFHandler ndofs=16
+DEAL:0::*** 1. CELLS ***
+DEAL:0::*** 2. FACES ***
+DEAL:0::* own_cells=0 ghost_cells=0 own_faces=1 faces_to_ghost=0
+DEAL:0::F cell1 = 0_2:00 face = 1 cell2 = 0_2:01 face2 = 0
+DEAL:0::F cell1 = 0_2:00 face = 3 cell2 = 0_2:02 face2 = 2
+DEAL:0::F cell1 = 0_2:01 face = 3 cell2 = 0_2:03 face2 = 2
+DEAL:0::F cell1 = 0_2:02 face = 1 cell2 = 0_2:03 face2 = 0
+DEAL:0::* own_cells=0 ghost_cells=0 own_faces=2 faces_to_ghost=0
+DEAL:0::F cell1 = 0_2:00 face = 1 cell2 = 0_2:01 face2 = 0
+DEAL:0::F cell1 = 0_2:00 face = 3 cell2 = 0_2:02 face2 = 2
+DEAL:0::F cell1 = 0_2:01 face = 0 cell2 = 0_2:00 face2 = 1
+DEAL:0::F cell1 = 0_2:01 face = 3 cell2 = 0_2:03 face2 = 2
+DEAL:0::F cell1 = 0_2:02 face = 1 cell2 = 0_2:03 face2 = 0
+DEAL:0::F cell1 = 0_2:02 face = 2 cell2 = 0_2:00 face2 = 3
+DEAL:0::F cell1 = 0_2:03 face = 0 cell2 = 0_2:02 face2 = 1
+DEAL:0::F cell1 = 0_2:03 face = 2 cell2 = 0_2:01 face2 = 3
+DEAL:0::* own_cells=0 ghost_cells=0 own_faces=0 faces_to_ghost=1
+DEAL:0::F cell1 = 0_2:01 face = 1 cell2 = 0_2:10 face2 = 0
+DEAL:0::F cell1 = 0_2:02 face = 3 cell2 = 0_2:20 face2 = 2
+DEAL:0::F cell1 = 0_2:03 face = 1 cell2 = 0_2:12 face2 = 0
+DEAL:0::F cell1 = 0_2:03 face = 3 cell2 = 0_2:21 face2 = 2
+DEAL:0::* own_cells=0 ghost_cells=0 own_faces=0 faces_to_ghost=2
+DEAL:0::F cell1 = 0_2:01 face = 1 cell2 = 0_2:10 face2 = 0
+DEAL:0::F cell1 = 0_2:02 face = 3 cell2 = 0_2:20 face2 = 2
+DEAL:0::F cell1 = 0_2:03 face = 1 cell2 = 0_2:12 face2 = 0
+DEAL:0::F cell1 = 0_2:03 face = 3 cell2 = 0_2:21 face2 = 2
+
+DEAL:1::DoFHandler ndofs=16
+DEAL:1::*** 1. CELLS ***
+DEAL:1::*** 2. FACES ***
+DEAL:1::* own_cells=0 ghost_cells=0 own_faces=1 faces_to_ghost=0
+DEAL:1::F cell1 = 0_2:10 face = 1 cell2 = 0_2:11 face2 = 0
+DEAL:1::F cell1 = 0_2:10 face = 3 cell2 = 0_2:12 face2 = 2
+DEAL:1::F cell1 = 0_2:11 face = 3 cell2 = 0_2:13 face2 = 2
+DEAL:1::F cell1 = 0_2:12 face = 1 cell2 = 0_2:13 face2 = 0
+DEAL:1::F cell1 = 0_2:20 face = 1 cell2 = 0_2:21 face2 = 0
+DEAL:1::F cell1 = 0_2:20 face = 3 cell2 = 0_2:22 face2 = 2
+DEAL:1::F cell1 = 0_2:21 face = 3 cell2 = 0_2:23 face2 = 2
+DEAL:1::F cell1 = 0_2:22 face = 1 cell2 = 0_2:23 face2 = 0
+DEAL:1::* own_cells=0 ghost_cells=0 own_faces=2 faces_to_ghost=0
+DEAL:1::F cell1 = 0_2:10 face = 1 cell2 = 0_2:11 face2 = 0
+DEAL:1::F cell1 = 0_2:10 face = 3 cell2 = 0_2:12 face2 = 2
+DEAL:1::F cell1 = 0_2:11 face = 0 cell2 = 0_2:10 face2 = 1
+DEAL:1::F cell1 = 0_2:11 face = 3 cell2 = 0_2:13 face2 = 2
+DEAL:1::F cell1 = 0_2:12 face = 1 cell2 = 0_2:13 face2 = 0
+DEAL:1::F cell1 = 0_2:12 face = 2 cell2 = 0_2:10 face2 = 3
+DEAL:1::F cell1 = 0_2:13 face = 0 cell2 = 0_2:12 face2 = 1
+DEAL:1::F cell1 = 0_2:13 face = 2 cell2 = 0_2:11 face2 = 3
+DEAL:1::F cell1 = 0_2:20 face = 1 cell2 = 0_2:21 face2 = 0
+DEAL:1::F cell1 = 0_2:20 face = 3 cell2 = 0_2:22 face2 = 2
+DEAL:1::F cell1 = 0_2:21 face = 0 cell2 = 0_2:20 face2 = 1
+DEAL:1::F cell1 = 0_2:21 face = 3 cell2 = 0_2:23 face2 = 2
+DEAL:1::F cell1 = 0_2:22 face = 1 cell2 = 0_2:23 face2 = 0
+DEAL:1::F cell1 = 0_2:22 face = 2 cell2 = 0_2:20 face2 = 3
+DEAL:1::F cell1 = 0_2:23 face = 0 cell2 = 0_2:22 face2 = 1
+DEAL:1::F cell1 = 0_2:23 face = 2 cell2 = 0_2:21 face2 = 3
+DEAL:1::* own_cells=0 ghost_cells=0 own_faces=0 faces_to_ghost=1
+DEAL:1::F cell1 = 0_2:12 face = 3 cell2 = 0_2:30 face2 = 2
+DEAL:1::F cell1 = 0_2:13 face = 3 cell2 = 0_2:31 face2 = 2
+DEAL:1::F cell1 = 0_2:21 face = 1 cell2 = 0_2:30 face2 = 0
+DEAL:1::F cell1 = 0_2:23 face = 1 cell2 = 0_2:32 face2 = 0
+DEAL:1::* own_cells=0 ghost_cells=0 own_faces=0 faces_to_ghost=2
+DEAL:1::F cell1 = 0_2:10 face = 0 cell2 = 0_2:01 face2 = 1
+DEAL:1::F cell1 = 0_2:12 face = 0 cell2 = 0_2:03 face2 = 1
+DEAL:1::F cell1 = 0_2:12 face = 3 cell2 = 0_2:30 face2 = 2
+DEAL:1::F cell1 = 0_2:13 face = 3 cell2 = 0_2:31 face2 = 2
+DEAL:1::F cell1 = 0_2:20 face = 2 cell2 = 0_2:02 face2 = 3
+DEAL:1::F cell1 = 0_2:21 face = 1 cell2 = 0_2:30 face2 = 0
+DEAL:1::F cell1 = 0_2:21 face = 2 cell2 = 0_2:03 face2 = 3
+DEAL:1::F cell1 = 0_2:23 face = 1 cell2 = 0_2:32 face2 = 0
+
+
+DEAL:2::DoFHandler ndofs=16
+DEAL:2::*** 1. CELLS ***
+DEAL:2::*** 2. FACES ***
+DEAL:2::* own_cells=0 ghost_cells=0 own_faces=1 faces_to_ghost=0
+DEAL:2::F cell1 = 0_2:30 face = 1 cell2 = 0_2:31 face2 = 0
+DEAL:2::F cell1 = 0_2:30 face = 3 cell2 = 0_2:32 face2 = 2
+DEAL:2::F cell1 = 0_2:31 face = 3 cell2 = 0_2:33 face2 = 2
+DEAL:2::F cell1 = 0_2:32 face = 1 cell2 = 0_2:33 face2 = 0
+DEAL:2::* own_cells=0 ghost_cells=0 own_faces=2 faces_to_ghost=0
+DEAL:2::F cell1 = 0_2:30 face = 1 cell2 = 0_2:31 face2 = 0
+DEAL:2::F cell1 = 0_2:30 face = 3 cell2 = 0_2:32 face2 = 2
+DEAL:2::F cell1 = 0_2:31 face = 0 cell2 = 0_2:30 face2 = 1
+DEAL:2::F cell1 = 0_2:31 face = 3 cell2 = 0_2:33 face2 = 2
+DEAL:2::F cell1 = 0_2:32 face = 1 cell2 = 0_2:33 face2 = 0
+DEAL:2::F cell1 = 0_2:32 face = 2 cell2 = 0_2:30 face2 = 3
+DEAL:2::F cell1 = 0_2:33 face = 0 cell2 = 0_2:32 face2 = 1
+DEAL:2::F cell1 = 0_2:33 face = 2 cell2 = 0_2:31 face2 = 3
+DEAL:2::* own_cells=0 ghost_cells=0 own_faces=0 faces_to_ghost=1
+DEAL:2::* own_cells=0 ghost_cells=0 own_faces=0 faces_to_ghost=2
+DEAL:2::F cell1 = 0_2:30 face = 0 cell2 = 0_2:21 face2 = 1
+DEAL:2::F cell1 = 0_2:30 face = 2 cell2 = 0_2:12 face2 = 3
+DEAL:2::F cell1 = 0_2:31 face = 2 cell2 = 0_2:13 face2 = 3
+DEAL:2::F cell1 = 0_2:32 face = 0 cell2 = 0_2:23 face2 = 1
+