Reference documentation for deal.II version 8.5.1
shared_tria.h
1 // ---------------------------------------------------------------------
2 // @f$Id: tria.h 32739 2014-04-08 16:39:47Z denis.davydov @f$
3 //
4 // Copyright (C) 2008 - 2016 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 #ifndef dealii__distributed__shared_tria_h
18 #define dealii__distributed__shared_tria_h
19 
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/subscriptor.h>
23 #include <deal.II/base/smartpointer.h>
24 #include <deal.II/base/template_constraints.h>
25 #include <deal.II/grid/tria.h>
26 
27 #include <deal.II/distributed/tria_base.h>
28 
29 #include <deal.II/base/std_cxx1x/function.h>
30 #include <deal.II/base/std_cxx1x/tuple.h>
31 
32 #include <set>
33 #include <vector>
34 #include <list>
35 #include <utility>
36 
37 #ifdef DEAL_II_WITH_MPI
38 # include <mpi.h>
39 #endif
40 
41 
42 DEAL_II_NAMESPACE_OPEN
43 
44 namespace parallel
45 {
46 
47 #ifdef DEAL_II_WITH_MPI
48 
49 
50  namespace shared
51  {
52 
67  template <int dim, int spacedim = dim>
68  class Triangulation : public ::parallel::Triangulation<dim,spacedim>
69  {
70  public:
71  typedef typename ::Triangulation<dim,spacedim>::active_cell_iterator active_cell_iterator;
72  typedef typename ::Triangulation<dim,spacedim>::cell_iterator cell_iterator;
73 
84  const typename ::Triangulation<dim,spacedim>::MeshSmoothing =
86  const bool allow_artificial_cells = false);
87 
91  virtual ~Triangulation ();
92 
101  virtual void execute_coarsening_and_refinement ();
102 
109  virtual void create_triangulation (const std::vector< Point< spacedim > > &vertices,
110  const std::vector< CellData< dim > > &cells,
111  const SubCellData &subcelldata);
112 
123  virtual void copy_triangulation (const ::Triangulation<dim, spacedim> &other_tria);
124 
133  template <class Archive>
134  void load (Archive &ar, const unsigned int version);
135 
144  const std::vector<types::subdomain_id> &get_true_subdomain_ids_of_cells() const;
145 
150  bool with_artificial_cells() const;
151 
152  private:
157 
162  void partition();
163 
175  std::vector<types::subdomain_id> true_subdomain_ids_of_cells;
176  };
177 
178  template <int dim, int spacedim>
179  template <class Archive>
180  void
182  const unsigned int version)
183  {
185  partition();
186  this->update_number_cache ();
187  }
188  }
189 #else
190 
191  namespace shared
192  {
193 
204  template <int dim, int spacedim = dim>
205  class Triangulation : public ::parallel::Triangulation<dim,spacedim>
206  {
207  public:
208 
212  const std::vector<types::subdomain_id> &get_true_subdomain_ids_of_cells() const;
213 
217  bool with_artificial_cells() const;
218  private:
222  Triangulation ();
223 
227  std::vector<types::subdomain_id> true_subdomain_ids_of_cells;
228  };
229  }
230 
231 
232 #endif
233 }
234 
235 DEAL_II_NAMESPACE_CLOSE
236 
237 #endif
virtual void execute_coarsening_and_refinement()
Definition: shared_tria.cc:117
Triangulation(MPI_Comm mpi_communicator, const typename ::Triangulation< dim, spacedim >::MeshSmoothing=(::Triangulation< dim, spacedim >::none), const bool allow_artificial_cells=false)
Definition: shared_tria.cc:37
std::vector< Point< spacedim > > vertices
Definition: tria.h:3447
const std::vector< types::subdomain_id > & get_true_subdomain_ids_of_cells() const
Definition: shared_tria.cc:102
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
Definition: shared_tria.cc:128
virtual void copy_triangulation(const ::Triangulation< dim, spacedim > &other_tria)
Definition: shared_tria.cc:152
void load(Archive &ar, const unsigned int version)
std::vector< types::subdomain_id > true_subdomain_ids_of_cells
Definition: shared_tria.h:175
void load(Archive &ar, const unsigned int version)
Definition: shared_tria.h:181