// $Id$
// Version: $Name$
//
-// Copyright (C) 2001, 2002, 2003, 2004, 2005 by the deal.II authors
+// Copyright (C) 2001, 2002, 2003, 2004, 2005, 2006 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
#include <grid/tria_accessor.h>
#include <grid/tria_iterator.h>
+#include <list>
* type of the first parameter
* may be either
* Triangulation,
- * DoFHandler, or
+ * DoFHandler, hp::DoFHandler, or
* MGDoFHandler, i.e. we
* can find the cell around a
* point for iterators into each
count_cells_with_subdomain_association (const Triangulation<dim> &triangulation,
const unsigned int subdomain);
+ /**
+ * Given two mesh containers
+ * (i.e. objects of type
+ * Triangulation, DoFHandler,
+ * hp::DoFHandler, or
+ * MGDoFHandler) that are based
+ * on the same coarse mesh, this
+ * function figures out a set of
+ * cells that are matched between
+ * the two meshes and where at
+ * most one of the meshes is more
+ * refined on this cell. In other
+ * words, it finds the smallest
+ * cells that are common to both
+ * meshes, and that together
+ * completely cover the domain.
+ *
+ * This function is useful, for
+ * example, in time-dependent or
+ * nonlinear application, where
+ * one has to integrate a
+ * solution defined on one mesh
+ * (e.g., the one from the
+ * previous time step or
+ * nonlinear iteration) against
+ * the shape functions of another
+ * mesh (the next time step, the
+ * next nonlinear iteration). If,
+ * for example, the new mesh is
+ * finer, then one has to obtain
+ * the solution on the coarse
+ * mesh (mesh_1) and interpolate
+ * it to the children of the
+ * corresponding cell of
+ * mesh_2. Conversely, if the new
+ * mesh is coarser, one has to
+ * express the coarse cell shape
+ * function by a linear
+ * combination of fine cell shape
+ * functions. In either case, one
+ * needs to loop over the finest
+ * cells that are common to both
+ * triangulations. This function
+ * returns a list of pairs of
+ * matching iterators to cells in
+ * the two meshes that can be
+ * used to this end.
+ *
+ * Note that the list of these
+ * iterators is not necessarily
+ * order, and does also not
+ * necessarily coincide with the
+ * order in which cells are
+ * traversed in one, or both, of
+ * the meshes given as arguments.
+ */
+ template <typename Container>
+ static
+ std::list<std::pair<typename Container::cell_iterator,
+ typename Container::cell_iterator> >
+ get_finest_common_cells (const Container &mesh_1,
+ const Container &mesh_2);
+
+ /**
+ * Return true if the two
+ * triangulations are based on
+ * the same coarse mesh. This is
+ * determined by checking whether
+ * they have the same number of
+ * cells on the coarsest level,
+ * and then checking that they
+ * have the same vertices.
+ *
+ * The two meshes may have
+ * different refinement histories
+ * beyond the coarse mesh.
+ */
+ template <int dim>
+ static
+ bool
+ have_same_coarse_mesh (const Triangulation<dim> &mesh_1,
+ const Triangulation<dim> &mesh_2);
+
+ /**
+ * The same function as above,
+ * but working on arguments of
+ * type DoFHandler,
+ * hp::DoFHandler, or
+ * MGDoFHandler. This function is
+ * provided to allow calling
+ * have_same_coarse_mesh for all
+ * types of containers
+ * representing triangulations or
+ * the classes built on
+ * triangulations.
+ */
+ template <typename Container>
+ static
+ bool
+ have_same_coarse_mesh (const Container &mesh_1,
+ const Container &mesh_2);
+
/**
* Exception
*/
#include <cmath>
+
#if deal_II_dimension != 1
template <int dim>
+template <typename Container>
+std::list<std::pair<typename Container::cell_iterator,
+ typename Container::cell_iterator> >
+GridTools::get_finest_common_cells (const Container &mesh_1,
+ const Container &mesh_2)
+{
+ Assert (have_same_coarse_mesh (mesh_1, mesh_2),
+ ExcMessage ("The two containers must be represent triangulations that "
+ "have the same coarse meshes"));
+
+ const unsigned int dim = Container::dimension;
+
+ // the algorithm goes as follows:
+ // first, we fill a list with pairs
+ // of iterators common to the two
+ // meshes on the coarsest
+ // level. then we traverse the
+ // list; each time, we find a pair
+ // of iterators for which both
+ // correspond to non-active cells,
+ // we delete this item and push the
+ // pairs of iterators to their
+ // children to the back. if these
+ // again both correspond to
+ // non-active cells, we will get to
+ // the later on for further
+ // consideration
+ typedef
+ std::list<std::pair<typename Container::cell_iterator,
+ typename Container::cell_iterator> >
+ CellList;
+
+ CellList cell_list;
+
+ // first push the coarse level cells
+ typename Container::cell_iterator
+ cell_1 = mesh_1.begin(0),
+ cell_2 = mesh_2.begin(0);
+ for (; cell_1 != mesh_1.end(0); ++cell_1, ++cell_2)
+ cell_list.push_back (std::make_pair (cell_1, cell_2));
+
+ // then traverse list as described
+ // above
+ typename CellList::iterator cell_pair = cell_list.begin();
+ while (cell_pair != cell_list.end())
+ {
+ // if both cells in this pair
+ // have children, then erase
+ // this element and push their
+ // children instead
+ if (cell_pair->first->has_children()
+ &&
+ cell_pair->second->has_children())
+ {
+ for (unsigned int c=0; c<GeometryInfo<dim>::children_per_cell; ++c)
+ cell_list.push_back (std::make_pair (cell_pair->first->child(c),
+ cell_pair->second->child(c)));
+
+ // erasing an iterator
+ // keeps other iterators
+ // valid, so already
+ // advance the present
+ // iterator by one and then
+ // delete the element we've
+ // visited before
+ const typename CellList::iterator previous_cell_pair = cell_pair;
+ ++cell_pair;
+
+ cell_list.erase (previous_cell_pair);
+ }
+ else
+ // both cells are active, do
+ // nothing
+ ++cell_pair;
+ }
+
+ // just to make sure everything is
+ // ok, validate that all pairs have
+ // at least one active iterator
+ for (cell_pair = cell_list.begin(); cell_pair != cell_list.end(); ++cell_pair)
+ Assert (!cell_pair->first->has_children()
+ ||
+ !cell_pair->second->has_children(),
+ ExcInternalError());
+
+ return cell_list;
+}
+
+
+
+template <int dim>
+bool
+GridTools::have_same_coarse_mesh (const Triangulation<dim> &mesh_1,
+ const Triangulation<dim> &mesh_2)
+{
+ // make sure the two meshes have
+ // the same number of coarse cells
+ if (mesh_1.n_cells (0) != mesh_2.n_cells (0))
+ return false;
+
+ // if so, also make sure they have
+ // the same vertices on the cells
+ // of the coarse mesh
+ typename Triangulation<dim>::cell_iterator
+ cell_1 = mesh_1.begin(0),
+ cell_2 = mesh_2.begin(0),
+ endc = mesh_1.end(0);
+ for (; cell_1!=endc; ++cell_1, ++cell_2)
+ for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
+ if (cell_1->vertex(v) != cell_2->vertex(v))
+ return false;
+
+ // if we've gotten through all
+ // this, then the meshes really
+ // seem to have a common coarse
+ // mesh
+ return true;
+}
+
+
+
+template <typename Container>
+bool
+GridTools::have_same_coarse_mesh (const Container &mesh_1,
+ const Container &mesh_2)
+{
+ return have_same_coarse_mesh (mesh_1.get_tria(),
+ mesh_2.get_tria());
+}
+
+
+
+
+// explicit instantiations
+
#if deal_II_dimension != 1
template
double
count_cells_with_subdomain_association (const Triangulation<deal_II_dimension> &,
const unsigned int );
+
+template
+std::list<std::pair<Triangulation<deal_II_dimension>::cell_iterator,
+ Triangulation<deal_II_dimension>::cell_iterator> >
+GridTools::
+get_finest_common_cells (const Triangulation<deal_II_dimension> &mesh_1,
+ const Triangulation<deal_II_dimension> &mesh_2);
+
+template
+std::list<std::pair<DoFHandler<deal_II_dimension>::cell_iterator,
+ DoFHandler<deal_II_dimension>::cell_iterator> >
+GridTools::
+get_finest_common_cells (const DoFHandler<deal_II_dimension> &mesh_1,
+ const DoFHandler<deal_II_dimension> &mesh_2);
+
+template
+std::list<std::pair<hp::DoFHandler<deal_II_dimension>::cell_iterator,
+ hp::DoFHandler<deal_II_dimension>::cell_iterator> >
+GridTools::
+get_finest_common_cells (const hp::DoFHandler<deal_II_dimension> &mesh_1,
+ const hp::DoFHandler<deal_II_dimension> &mesh_2);
+
+template
+std::list<std::pair<MGDoFHandler<deal_II_dimension>::cell_iterator,
+ MGDoFHandler<deal_II_dimension>::cell_iterator> >
+GridTools::
+get_finest_common_cells (const MGDoFHandler<deal_II_dimension> &mesh_1,
+ const MGDoFHandler<deal_II_dimension> &mesh_2);
+
+
+template
+bool
+GridTools::have_same_coarse_mesh (const Triangulation<deal_II_dimension> &mesh_1,
+ const Triangulation<deal_II_dimension> &mesh_2);
+
+template
+bool
+GridTools::have_same_coarse_mesh (const DoFHandler<deal_II_dimension> &mesh_1,
+ const DoFHandler<deal_II_dimension> &mesh_2);
+
+template
+bool
+GridTools::have_same_coarse_mesh (const hp::DoFHandler<deal_II_dimension> &mesh_1,
+ const hp::DoFHandler<deal_II_dimension> &mesh_2);
+
+template
+bool
+GridTools::have_same_coarse_mesh (const MGDoFHandler<deal_II_dimension> &mesh_1,
+ const MGDoFHandler<deal_II_dimension> &mesh_2);
<h3>deal.II</h3>
<ol>
+ <li> <p>
+ New: The function <code>GridTools::get_finest_common_cells</code> can be
+ used to find those cells that are shared by two triangulations that are
+ based on the same coarse mesh. This is useful, for example, when having
+ to integrate the solution on one mesh against the shape functions on
+ another mesh, a situation that frequently happens in time-dependent but
+ also in nonlinear problems.
+ <br>
+ (WB 2006/04/14)
+ </p>
+
+ <li> <p>
+ New: The function <code>GridTools::have_same_coarse_mesh</code> figures
+ out whether two triangulations, two DoFHandlers, etc, are built on the
+ same coarse mesh.
+ <br>
+ (WB 2006/04/14)
+ </p>
+
<li> <p>
New: Since calling <code>cell->get_dof_indices</code> is a
fairly frequent operation (called about 7 times per cell in
face_orientations_3d \
mg_dof_handler \
subcelldata \
- interpolate_boundary_values_*
+ interpolate_boundary_values_* \
+ have_same_coarse_mesh_* \
+ get_finest_common_cells_*
# from above list of regular expressions, generate the real set of
# tests
--- /dev/null
+//---------------------------- get_finest_common_cells_01.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- get_finest_common_cells_01.cc ---------------------------
+// check GridTools::get_finest_common_cells for Triangulation arguments
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+#include <grid/tria_accessor.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 2 triangulations with the
+ // same coarse grid, and refine
+ // them differently
+ Triangulation<dim> tria[2];
+
+ GridGenerator::hyper_cube (tria[0]);
+ GridGenerator::hyper_cube (tria[1]);
+
+ tria[0].refine_global (2);
+ tria[1].refine_global (2);
+
+ tria[0].begin_active()->set_refine_flag();
+ tria[0].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ typedef
+ std::list<std::pair<typename Triangulation<dim>::cell_iterator,
+ typename Triangulation<dim>::cell_iterator> >
+ CellList;
+
+ const CellList cell_list
+ = GridTools::get_finest_common_cells (tria[0], tria[1]);
+ for (typename CellList::const_iterator cell_pair = cell_list.begin();
+ cell_pair != cell_list.end(); ++cell_pair)
+ deallog << cell_pair->first << ' ' << cell_pair->second
+ << std::endl;
+}
+
+
+int main()
+{
+ std::ofstream logfile ("get_finest_common_cells_01/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.16 2.16
+DEAL::2.17 2.17
+DEAL::2.18 2.18
+DEAL::2.19 2.19
+DEAL::2.20 2.20
+DEAL::2.21 2.21
+DEAL::2.22 2.22
+DEAL::2.23 2.23
+DEAL::2.24 2.24
+DEAL::2.25 2.25
+DEAL::2.26 2.26
+DEAL::2.27 2.27
+DEAL::2.28 2.28
+DEAL::2.29 2.29
+DEAL::2.30 2.30
+DEAL::2.31 2.31
+DEAL::2.32 2.32
+DEAL::2.33 2.33
+DEAL::2.34 2.34
+DEAL::2.35 2.35
+DEAL::2.36 2.36
+DEAL::2.37 2.37
+DEAL::2.38 2.38
+DEAL::2.39 2.39
+DEAL::2.40 2.40
+DEAL::2.41 2.41
+DEAL::2.42 2.42
+DEAL::2.43 2.43
+DEAL::2.44 2.44
+DEAL::2.45 2.45
+DEAL::2.46 2.46
+DEAL::2.47 2.47
+DEAL::2.48 2.48
+DEAL::2.49 2.49
+DEAL::2.50 2.50
+DEAL::2.51 2.51
+DEAL::2.52 2.52
+DEAL::2.53 2.53
+DEAL::2.54 2.54
+DEAL::2.55 2.55
+DEAL::2.56 2.56
+DEAL::2.57 2.57
+DEAL::2.58 2.58
+DEAL::2.59 2.59
+DEAL::2.60 2.60
+DEAL::2.61 2.61
+DEAL::2.62 2.62
+DEAL::2.63 2.63
--- /dev/null
+//---------------------------- get_finest_common_cells_02.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- get_finest_common_cells_02.cc ---------------------------
+// check GridTools::get_finest_common_cells for DoFHandler arguments
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+#include <dofs/dof_accessor.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 2 triangulations with the
+ // same coarse grid, and refine
+ // them differently
+ Triangulation<dim> tria[2];
+
+ GridGenerator::hyper_cube (tria[0]);
+ GridGenerator::hyper_cube (tria[1]);
+
+ tria[0].refine_global (2);
+ tria[1].refine_global (2);
+
+ tria[0].begin_active()->set_refine_flag();
+ tria[0].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ DoFHandler<dim> dh0 (tria[0]);
+ DoFHandler<dim> dh1 (tria[1]);
+
+ typedef
+ std::list<std::pair<typename DoFHandler<dim>::cell_iterator,
+ typename DoFHandler<dim>::cell_iterator> >
+ CellList;
+
+ const CellList cell_list
+ = GridTools::get_finest_common_cells (dh0, dh1);
+ for (typename CellList::const_iterator cell_pair = cell_list.begin();
+ cell_pair != cell_list.end(); ++cell_pair)
+ deallog << cell_pair->first << ' ' << cell_pair->second
+ << std::endl;
+}
+
+
+int main()
+{
+ std::ofstream logfile ("get_finest_common_cells_02/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.16 2.16
+DEAL::2.17 2.17
+DEAL::2.18 2.18
+DEAL::2.19 2.19
+DEAL::2.20 2.20
+DEAL::2.21 2.21
+DEAL::2.22 2.22
+DEAL::2.23 2.23
+DEAL::2.24 2.24
+DEAL::2.25 2.25
+DEAL::2.26 2.26
+DEAL::2.27 2.27
+DEAL::2.28 2.28
+DEAL::2.29 2.29
+DEAL::2.30 2.30
+DEAL::2.31 2.31
+DEAL::2.32 2.32
+DEAL::2.33 2.33
+DEAL::2.34 2.34
+DEAL::2.35 2.35
+DEAL::2.36 2.36
+DEAL::2.37 2.37
+DEAL::2.38 2.38
+DEAL::2.39 2.39
+DEAL::2.40 2.40
+DEAL::2.41 2.41
+DEAL::2.42 2.42
+DEAL::2.43 2.43
+DEAL::2.44 2.44
+DEAL::2.45 2.45
+DEAL::2.46 2.46
+DEAL::2.47 2.47
+DEAL::2.48 2.48
+DEAL::2.49 2.49
+DEAL::2.50 2.50
+DEAL::2.51 2.51
+DEAL::2.52 2.52
+DEAL::2.53 2.53
+DEAL::2.54 2.54
+DEAL::2.55 2.55
+DEAL::2.56 2.56
+DEAL::2.57 2.57
+DEAL::2.58 2.58
+DEAL::2.59 2.59
+DEAL::2.60 2.60
+DEAL::2.61 2.61
+DEAL::2.62 2.62
+DEAL::2.63 2.63
--- /dev/null
+//---------------------------- get_finest_common_cells_03.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- get_finest_common_cells_03.cc ---------------------------
+// check GridTools::get_finest_common_cells for hp::DoFHandler arguments
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+#include <dofs/dof_accessor.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 2 triangulations with the
+ // same coarse grid, and refine
+ // them differently
+ Triangulation<dim> tria[2];
+
+ GridGenerator::hyper_cube (tria[0]);
+ GridGenerator::hyper_cube (tria[1]);
+
+ tria[0].refine_global (2);
+ tria[1].refine_global (2);
+
+ tria[0].begin_active()->set_refine_flag();
+ tria[0].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ hp::DoFHandler<dim> dh0 (tria[0]);
+ hp::DoFHandler<dim> dh1 (tria[1]);
+
+ typedef
+ std::list<std::pair<typename hp::DoFHandler<dim>::cell_iterator,
+ typename hp::DoFHandler<dim>::cell_iterator> >
+ CellList;
+
+ const CellList cell_list
+ = GridTools::get_finest_common_cells (dh0, dh1);
+ for (typename CellList::const_iterator cell_pair = cell_list.begin();
+ cell_pair != cell_list.end(); ++cell_pair)
+ deallog << cell_pair->first << ' ' << cell_pair->second
+ << std::endl;
+}
+
+
+int main()
+{
+ std::ofstream logfile ("get_finest_common_cells_03/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.16 2.16
+DEAL::2.17 2.17
+DEAL::2.18 2.18
+DEAL::2.19 2.19
+DEAL::2.20 2.20
+DEAL::2.21 2.21
+DEAL::2.22 2.22
+DEAL::2.23 2.23
+DEAL::2.24 2.24
+DEAL::2.25 2.25
+DEAL::2.26 2.26
+DEAL::2.27 2.27
+DEAL::2.28 2.28
+DEAL::2.29 2.29
+DEAL::2.30 2.30
+DEAL::2.31 2.31
+DEAL::2.32 2.32
+DEAL::2.33 2.33
+DEAL::2.34 2.34
+DEAL::2.35 2.35
+DEAL::2.36 2.36
+DEAL::2.37 2.37
+DEAL::2.38 2.38
+DEAL::2.39 2.39
+DEAL::2.40 2.40
+DEAL::2.41 2.41
+DEAL::2.42 2.42
+DEAL::2.43 2.43
+DEAL::2.44 2.44
+DEAL::2.45 2.45
+DEAL::2.46 2.46
+DEAL::2.47 2.47
+DEAL::2.48 2.48
+DEAL::2.49 2.49
+DEAL::2.50 2.50
+DEAL::2.51 2.51
+DEAL::2.52 2.52
+DEAL::2.53 2.53
+DEAL::2.54 2.54
+DEAL::2.55 2.55
+DEAL::2.56 2.56
+DEAL::2.57 2.57
+DEAL::2.58 2.58
+DEAL::2.59 2.59
+DEAL::2.60 2.60
+DEAL::2.61 2.61
+DEAL::2.62 2.62
+DEAL::2.63 2.63
--- /dev/null
+//---------------------------- get_finest_common_cells_04.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- get_finest_common_cells_04.cc ---------------------------
+// check GridTools::get_finest_common_cells for MGDoFHandler arguments
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+#include <multigrid/mg_dof_accessor.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 2 triangulations with the
+ // same coarse grid, and refine
+ // them differently
+ Triangulation<dim> tria[2];
+
+ GridGenerator::hyper_cube (tria[0]);
+ GridGenerator::hyper_cube (tria[1]);
+
+ tria[0].refine_global (2);
+ tria[1].refine_global (2);
+
+ tria[0].begin_active()->set_refine_flag();
+ tria[0].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ tria[1].last_active()->set_refine_flag();
+ tria[1].execute_coarsening_and_refinement ();
+
+ MGDoFHandler<dim> dh0 (tria[0]);
+ MGDoFHandler<dim> dh1 (tria[1]);
+
+ typedef
+ std::list<std::pair<typename MGDoFHandler<dim>::cell_iterator,
+ typename MGDoFHandler<dim>::cell_iterator> >
+ CellList;
+
+ const CellList cell_list
+ = GridTools::get_finest_common_cells (dh0, dh1);
+ for (typename CellList::const_iterator cell_pair = cell_list.begin();
+ cell_pair != cell_list.end(); ++cell_pair)
+ deallog << cell_pair->first << ' ' << cell_pair->second
+ << std::endl;
+}
+
+
+int main()
+{
+ std::ofstream logfile ("get_finest_common_cells_04/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.0 2.0
+DEAL::2.1 2.1
+DEAL::2.2 2.2
+DEAL::2.3 2.3
+DEAL::2.4 2.4
+DEAL::2.5 2.5
+DEAL::2.6 2.6
+DEAL::2.7 2.7
+DEAL::2.8 2.8
+DEAL::2.9 2.9
+DEAL::2.10 2.10
+DEAL::2.11 2.11
+DEAL::2.12 2.12
+DEAL::2.13 2.13
+DEAL::2.14 2.14
+DEAL::2.15 2.15
+DEAL::2.16 2.16
+DEAL::2.17 2.17
+DEAL::2.18 2.18
+DEAL::2.19 2.19
+DEAL::2.20 2.20
+DEAL::2.21 2.21
+DEAL::2.22 2.22
+DEAL::2.23 2.23
+DEAL::2.24 2.24
+DEAL::2.25 2.25
+DEAL::2.26 2.26
+DEAL::2.27 2.27
+DEAL::2.28 2.28
+DEAL::2.29 2.29
+DEAL::2.30 2.30
+DEAL::2.31 2.31
+DEAL::2.32 2.32
+DEAL::2.33 2.33
+DEAL::2.34 2.34
+DEAL::2.35 2.35
+DEAL::2.36 2.36
+DEAL::2.37 2.37
+DEAL::2.38 2.38
+DEAL::2.39 2.39
+DEAL::2.40 2.40
+DEAL::2.41 2.41
+DEAL::2.42 2.42
+DEAL::2.43 2.43
+DEAL::2.44 2.44
+DEAL::2.45 2.45
+DEAL::2.46 2.46
+DEAL::2.47 2.47
+DEAL::2.48 2.48
+DEAL::2.49 2.49
+DEAL::2.50 2.50
+DEAL::2.51 2.51
+DEAL::2.52 2.52
+DEAL::2.53 2.53
+DEAL::2.54 2.54
+DEAL::2.55 2.55
+DEAL::2.56 2.56
+DEAL::2.57 2.57
+DEAL::2.58 2.58
+DEAL::2.59 2.59
+DEAL::2.60 2.60
+DEAL::2.61 2.61
+DEAL::2.62 2.62
+DEAL::2.63 2.63
--- /dev/null
+//---------------------------- have_same_coarse_mesh_01.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- have_same_coarse_mesh_01.cc ---------------------------
+// check GridTools::have_same_coarse_mesh for triangulations
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 3 triangulations
+ Triangulation<dim> tria[3];
+
+ GridGenerator::hyper_cube (tria[0]);
+ tria[0].refine_global (1);
+
+ GridGenerator::hyper_cube (tria[1]);
+ GridTools::scale (2, tria[1]);
+ tria[1].refine_global (2);
+
+ if (dim != 1)
+ GridGenerator::hyper_ball (tria[2]);
+ else
+ {
+ GridGenerator::hyper_cube (tria[2]);
+ GridTools::shift (Point<dim>(2.), tria[2]);
+ }
+
+ tria[2].refine_global (3);
+
+ for (unsigned int i=0; i<3; ++i)
+ for (unsigned int j=0; j<3; ++j)
+ {
+ Assert (GridTools::have_same_coarse_mesh (tria[i], tria[j])
+ ==
+ (i == j),
+ ExcInternalError());
+
+ deallog << "meshes " << i << " and " << j << ": "
+ << (GridTools::have_same_coarse_mesh (tria[i], tria[j])
+ ?
+ "true"
+ :
+ "false")
+ << std::endl;
+ }
+}
+
+
+int main()
+{
+ std::ofstream logfile ("have_same_coarse_mesh_01/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
--- /dev/null
+//---------------------------- have_same_coarse_mesh_02.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- have_same_coarse_mesh_02.cc ---------------------------
+// check GridTools::have_same_coarse_mesh for DoFHandler arguments
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 3 triangulations
+ Triangulation<dim> tria[3];
+
+ GridGenerator::hyper_cube (tria[0]);
+ tria[0].refine_global (1);
+
+ GridGenerator::hyper_cube (tria[1]);
+ GridTools::scale (2, tria[1]);
+ tria[1].refine_global (2);
+
+ if (dim != 1)
+ GridGenerator::hyper_ball (tria[2]);
+ else
+ {
+ GridGenerator::hyper_cube (tria[2]);
+ GridTools::shift (Point<dim>(2.), tria[2]);
+ }
+
+ tria[2].refine_global (3);
+
+ DoFHandler<dim> dh0 (tria[0]);
+ DoFHandler<dim> dh1 (tria[1]);
+ DoFHandler<dim> dh2 (tria[2]);
+
+ DoFHandler<dim> *dof_handler[3] = { &dh0, &dh1, &dh2 };
+
+ for (unsigned int i=0; i<3; ++i)
+ for (unsigned int j=0; j<3; ++j)
+ {
+ Assert (GridTools::have_same_coarse_mesh (*dof_handler[i], *dof_handler[j])
+ ==
+ (i == j),
+ ExcInternalError());
+
+ deallog << "meshes " << i << " and " << j << ": "
+ << (GridTools::have_same_coarse_mesh (*dof_handler[i], *dof_handler[j])
+ ?
+ "true"
+ :
+ "false")
+ << std::endl;
+ }
+}
+
+
+int main()
+{
+ std::ofstream logfile ("have_same_coarse_mesh_02/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
--- /dev/null
+//---------------------------- have_same_coarse_mesh_03.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- have_same_coarse_mesh_03.cc ---------------------------
+// check GridTools::have_same_coarse_mesh for hp::DoFHandler arguments
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 3 triangulations
+ Triangulation<dim> tria[3];
+
+ GridGenerator::hyper_cube (tria[0]);
+ tria[0].refine_global (1);
+
+ GridGenerator::hyper_cube (tria[1]);
+ GridTools::scale (2, tria[1]);
+ tria[1].refine_global (2);
+
+ if (dim != 1)
+ GridGenerator::hyper_ball (tria[2]);
+ else
+ {
+ GridGenerator::hyper_cube (tria[2]);
+ GridTools::shift (Point<dim>(2.), tria[2]);
+ }
+
+ tria[2].refine_global (3);
+
+ hp::DoFHandler<dim> dh0 (tria[0]);
+ hp::DoFHandler<dim> dh1 (tria[1]);
+ hp::DoFHandler<dim> dh2 (tria[2]);
+
+ hp::DoFHandler<dim> *dof_handler[3] = { &dh0, &dh1, &dh2 };
+
+ for (unsigned int i=0; i<3; ++i)
+ for (unsigned int j=0; j<3; ++j)
+ {
+ Assert (GridTools::have_same_coarse_mesh (*dof_handler[i], *dof_handler[j])
+ ==
+ (i == j),
+ ExcInternalError());
+
+ deallog << "meshes " << i << " and " << j << ": "
+ << (GridTools::have_same_coarse_mesh (*dof_handler[i], *dof_handler[j])
+ ?
+ "true"
+ :
+ "false")
+ << std::endl;
+ }
+}
+
+
+int main()
+{
+ std::ofstream logfile ("have_same_coarse_mesh_03/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
--- /dev/null
+//---------------------------- have_same_coarse_mesh_04.cc ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 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.
+//
+//---------------------------- have_same_coarse_mesh_04.cc ---------------------------
+// check GridTools::have_same_coarse_mesh for MGDoFHandler arguments
+
+
+#include "../tests.h"
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/grid_tools.h>
+#include <dofs/dof_handler.h>
+#include <dofs/hp_dof_handler.h>
+#include <multigrid/mg_dof_handler.h>
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+ // create 3 triangulations
+ Triangulation<dim> tria[3];
+
+ GridGenerator::hyper_cube (tria[0]);
+ tria[0].refine_global (1);
+
+ GridGenerator::hyper_cube (tria[1]);
+ GridTools::scale (2, tria[1]);
+ tria[1].refine_global (2);
+
+ if (dim != 1)
+ GridGenerator::hyper_ball (tria[2]);
+ else
+ {
+ GridGenerator::hyper_cube (tria[2]);
+ GridTools::shift (Point<dim>(2.), tria[2]);
+ }
+
+ tria[2].refine_global (3);
+
+ MGDoFHandler<dim> dh0 (tria[0]);
+ MGDoFHandler<dim> dh1 (tria[1]);
+ MGDoFHandler<dim> dh2 (tria[2]);
+
+ MGDoFHandler<dim> *dof_handler[3] = { &dh0, &dh1, &dh2 };
+
+ for (unsigned int i=0; i<3; ++i)
+ for (unsigned int j=0; j<3; ++j)
+ {
+ Assert (GridTools::have_same_coarse_mesh (*dof_handler[i], *dof_handler[j])
+ ==
+ (i == j),
+ ExcInternalError());
+
+ deallog << "meshes " << i << " and " << j << ": "
+ << (GridTools::have_same_coarse_mesh (*dof_handler[i], *dof_handler[j])
+ ?
+ "true"
+ :
+ "false")
+ << std::endl;
+ }
+}
+
+
+int main()
+{
+ std::ofstream logfile ("have_same_coarse_mesh_04/output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<1>();
+ test<2>();
+ test<3>();
+}
+
--- /dev/null
+
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true
+DEAL::meshes 0 and 0: true
+DEAL::meshes 0 and 1: false
+DEAL::meshes 0 and 2: false
+DEAL::meshes 1 and 0: false
+DEAL::meshes 1 and 1: true
+DEAL::meshes 1 and 2: false
+DEAL::meshes 2 and 0: false
+DEAL::meshes 2 and 1: false
+DEAL::meshes 2 and 2: true