Flow123d
local_to_global_map.hh
Go to the documentation of this file.
1 /*
2  * local_to_global_map.hh
3  *
4  * Created on: Mar 9, 2012
5  * Author: jb
6  */
7 
8 #ifndef LOCAL_TO_GLOBAL_MAP_HH_
9 #define LOCAL_TO_GLOBAL_MAP_HH_
10 
11 #include <set>
12 #include <vector>
13 
14 /// Using Boost shared pointer.
15 #include <boost/smart_ptr/shared_ptr.hpp>
16 #include <boost/smart_ptr/make_shared.hpp>
17 
18 #include "system/system.hh"
19 #include "system/global_defs.h"
20 
21 //namespace boost{
22 // template<class T> class shared_ptr;
23 //}
24 
25 class Distribution;
26 
27 /**
28  * @brief class to manage local indices on sub-domain to global indices on domain
29  *
30  * Currently (March 2012) the local to global maps are managed by individual equations (e.g. DarcyFlow has el_4_loc .. there it is complicated
31  * by different local indexing of elements, sides and edges). In fact el_4_loc etc. are new_local to old_global maps.
32  * This map is created as follows:
33  *
34  * PETSC solver:
35  * 1) create graph
36  * 2) make partitioning
37  * 2.5) create id_4_old (just used to create new_4_id, id was non continuous index not used anymore)
38  * 3) call id_maps which:
39  * 4) call ISPartitioninToNumbering : assign partitions to processors (identity, no optimization); make "distribution";
40  * make mapping (array of ints): old_local to new_global
41  * 5) from this we crate AO (application ordering from PETSc)
42  * 6) use AO to map identity array to old numbering -> creates map: new_global to old_global
43  * 7) go through new_local continuous part and map it to old_global index
44  * ---
45  * So this produce new_local to old_global mapping and only for continuous part of local indices (not for ghost indices)
46  * See that this works without any information about connectivity.
47  *
48  * METIS/BDDC solver:
49  * 1) graph, partition of elements, maps for elements (no overlap, no ghost values) using id_maps function
50  * 2) same for edges -> produce non overlapping local to global maps.
51  * 3) use std::set to collect all dofs on local elements (pass through mesh) complexity: n*log(n) n-local number of dofs
52  * 4) copy set to vector - makes local ordering arbitrary
53  *
54  * Future usage is:
55  * - in mesh to map local idx of entities to global indices
56  * - in dof handler to map local dofs idx to global
57  *
58  * In both cases the mapping is created by adding all global numbers on local subdomain. In order to keep local part of the map
59  * continuous, we need Distribution that describes splitting of global indices into continuous blocks. This Distribution is known
60  * as soon as we assign partitions to processors. So we assume that it is known at construction time.
61  *
62  */
63 
64 
65 
67 public:
68  /**
69  * Constructor. Starts filling of the map.
70  *
71  * @param distr Non overlapping distribution of global indices to processors in continuous blocks.
72  * Local block of indices forms first part of the mapping and then nonlocal indices follows.
73  * This constructor makes a deep copy of the distribution.
74  *
75  */
76  LocalToGlobalMap(const Distribution &distr);
77 
78  /**
79  * Same as the previous constructor, but just takes copy of shared pointer. This assumes that Distribution is allocated on the heap
80  * by something like:
81  *
82  * boost::smart_ptr<Distribution> distr(new Distribution(...));
83  */
84  LocalToGlobalMap(boost::shared_ptr<Distribution> distr);
85 
86  /**
87  * Insert a global index to the mapping.
88  */
89  void insert(const unsigned int idx);
90  /**
91  * Insert more indices at once.
92  */
93  void insert(const std::vector<unsigned int> &indices);
94  /**
95  * Finish filling stage. Creates mapping array, then you can use () operator to map indices.
96  */
97  void finalize();
98  /**
99  * Maps local index to the global one. For DEBUG, performs check for dimension.
100  */
101  inline unsigned int operator[] (const unsigned int local_idx) const
102  {
103  ASSERT_LESS( local_idx, global_indices_.size() );
104  return global_indices_[local_idx];
105  }
106 
107  /**
108  * Returns size of local map.
109  */
110  unsigned int size() const
111  { return global_indices_.size(); }
112 
113  /**
114  * Returns smart_ptr to the Distribution of the global indices. Allow share this among many objects.
115  */
116  boost::shared_ptr<Distribution> &get_distr()
117  { return distr_; }
118 
119  /**
120  * Returns inner vector.
121  */
123  { return global_indices_; }
124 
125 
126 
127 private:
128  /// auxiliary set of non-local part of the map
129  std::set<unsigned int> *nonlocal_indices_;
130  /// distribution of the global indices
131  boost::shared_ptr<Distribution> distr_;
132  /// mapping for all local indices
134 };
135 
136 
137 #endif /* LOCAL_TO_GLOBAL_MAP_HH_ */