Flow123d
region.hh
Go to the documentation of this file.
1 /*
2  * region.hh
3  *
4  * Created on: Nov 27, 2012
5  * Author: jb
6  *
7  * TODO:
8  * - komentar k RegionIdx
9  * - presun enum RegionType do public Region - komentar + pouzit v kodu
10  * - zkontrolovat chybove hlasky a ASSERTY, co z toho by melo byt pres exception?
11  *
12  * - Seems that GMSH allows repeating ID and Label on regions of different dimension, therefore
13  * label and ID are not unique without dimension.
14  *
15  */
16 
17 #ifndef REGION_HH_
18 #define REGION_HH_
19 
20 #include <string>
21 #include <vector>
22 #include <map>
23 
24 #include "system/system.hh"
25 #include "system/global_defs.h"
26 #include "system/exceptions.hh"
27 
28 #include "input/accessors.hh"
29 
30 #include <boost/multi_index_container.hpp>
31 #include <boost/multi_index/ordered_index.hpp>
32 #include <boost/multi_index/random_access_index.hpp>
33 #include <boost/multi_index/hashed_index.hpp>
34 #include <boost/multi_index/member.hpp>
35 
36 namespace BMI=::boost::multi_index;
37 
38 // forward declarations
39 class Region;
40 class RegionDB;
41 namespace Input {
42  namespace Type { class Record; }
43  class Record;
44  class Array;
45 }
46 
47 /**
48  * Base class that contains information about region:
49  * 1) contains integer value that specifies region
50  * 2) detects if region is valid
51  * 3) detects if the region is the bulk or boundary
52  * 4) detects boundary index of boundary region or bulk index of bulk region
53  */
54 class RegionIdx {
55 public:
56 
57  /// Default region is undefined/invalid
59 
60  /// Allow implicit conversion from Region. We loose information about input ID, label, dim stored in database.
61  //RegionIdx(Region region)
62  //: idx_(region.idx())
63  //{}
64 
65  /// Returns true if it is a Boundary region and false if it is a Bulk region.
66  inline bool is_boundary() const
67  { return !(idx_ & 1); }
68 
69  /// Returns false if the region has undefined/invalid value
70  inline bool is_valid() const
71  { return idx_!=undefined;}
72 
73  /// Returns a global index of the region.
74  inline unsigned int idx() const
75  { return idx_; }
76 
77  /// Returns index of the region in the boundary set.
78  inline unsigned int boundary_idx() const {
79  ASSERT( is_boundary(), "Try to get boundary index of a bulk region with internal idx: %d\n", idx_ );
80  return idx_ >> 1; }
81 
82  /// Returns index of the region in the bulk set.
83  inline unsigned int bulk_idx() const {
84  ASSERT( ! is_boundary(), "Try to get bulk index of boundary region with internal idx: %d\n", idx_ );
85  return idx_ >> 1; }
86 
87  /// Equality comparison operators for regions.
88  inline bool operator==(const RegionIdx &other) const
89  { return idx_ == other.idx_; }
90 
91  /// Equality comparison operators for regions.
92  inline bool operator!=(const RegionIdx &other) const
93  { return idx_ != other.idx_; }
94 
95 
96 protected:
97  /**
98  * Create accessor from the index. Should be private since implementation specific.
99  * We need some way how to iterate over: all regions, boundary regions, bulk regions -
100  * solution: have specific RegionSets for these three cases.
101  */
102  RegionIdx(unsigned int index)
103  : idx_(index) { }
104 
105  /**
106  * Internal region index. Regions of one RegionDB (corresponding to one mesh) forms more or less continuous sequence.
107  */
108  unsigned int idx_;
109 
110  /// index for undefined region
111  static const unsigned int undefined=0xffffffff;
112 
113  friend class Region;
114 };
115 
116 
117 
118 /**
119  * Class that represents disjoint part of computational domain (or domains). It consists of one integer value
120  * but provides access to other data stored in RegionDB. In particular provides string label and integer ID (unordered)
121  * further it provides fast (inlined) methods to:
122  * 1) detect if the region is the bulk region or boundary region
123  * 2) return index (this is used to select correct Field, possibly we can distinguish boundary_index and bulk_index)
124  *
125  * Implementation: currently we number bulk regions by odd indices and boundary regions by even indices.
126  *
127  */
128 class Region : public RegionIdx {
129 public:
130 
131  /**
132  * Types of region in mesh (bulk or boundary)
133  */
134  enum RegionType {
135  bulk=false,
136  boundary=true
137  };
138 
139 
140  /// Default region is undefined/invalid
142  : RegionIdx(), db_(NULL)
143  {}
144 
145  /// This should be used for construction from known RegionIdx. (e.g. in Mesh)
146  /// Do not use unless you can not get Region in other way.
147  Region(RegionIdx r_idx, const RegionDB & db)
148  : RegionIdx(r_idx), db_(&db)
149  {}
150 
152  {return RegionIdx(idx_); }
153 
154  /// Returns label of the region (using RegionDB)
155  std::string label() const;
156 
157  /// Returns id of the region (using RegionDB)
158  unsigned int id() const;
159 
160  /// Returns dimension of the region.
161  unsigned int dim() const;
162 
163  /* NOTE: seems that this conversion is provided automatically by C++, we taest this in region_test.cpp
164  * moreover the default version provided by compiler overrides our implementation
165  /// Conversion to RegionIdx class
166  /// NOTE:
167  //inline explicit operator RegionIdx() const {
168  // return RegionIdx(this->idx_);
169  //
170  */
171 
172  /**
173  * Returns static region database. Meant to be used for getting range of
174  * global, boundary, and bulk region indices.
175  */
176  const RegionDB &db() {
177  return *db_;
178  }
179 
180 private:
181  /**
182  * Create accessor from the index. Should be private since implementation specific.
183  * We need some way how to iterate over: all regions, boundary regions, bulk regions -
184  * solution: have specific RegionSets for these three cases.
185  */
186  Region(unsigned int index, const RegionDB &db)
187  : RegionIdx(index), db_(&db)
188  {}
189 
190  /// Comparative method of two regions
191  static bool comp(const Region &a, const Region &b)
192  { return a.idx_ < b.idx_; }
193 
194  /// Global variable with information about all regions.
195  const RegionDB *db_;
196 
197  friend class RegionDB;
198 };
199 
200 
201 
202 
203 /**
204  * Type representing a set of regions.
205  * CAn be used to set function(field) on more regions at once, possibly across meshes
206  *
207  * Regions stored in region set are always unique
208  */
210 
211 
212 
213 /**
214  * Class for conversion between an index and string label of an material.
215  * Class contains only static members in order to provide globally consistent
216  * indexing of materials across various meshes and functions.
217  *
218  * The conversion should be performed only through the input and output so
219  * that the lookup overhead could be shadowed by IO operations.
220  *
221  * Taking sizes and creating region sets should be possible only after the database is closed.
222  * We assume that all regions are known at beginning of the program (typically after reading all meshes)
223  * however they need not be used through the whole computation.
224  *
225  * TODO:
226  * In order to support more meshes , possibly changing during the time we need better policy for RegionDB closing.
227  * Currently we close RegionDB at first call to any of @p size methods. We need size information for initialization of
228  * RegionFields. However, every RegionField should be used for assembly over just one mesh (or set of meshes - supermesh?),
229  * surly this mesh has to be initialized before assembly so it could be initialized before initialization of RegionField which
230  * lives on this mesh. So the solution can be: RegionDB keeps list of meshes that has their regions registered in RegionDB.
231  * RegionField has signature of its mesh and check if the mesh is registered in the RegionDB before initialization of REgionField.
232  *
233  * Update TODO:
234  *
235  * - make one instance of Region DB per mesh; put region indices into mesh elements make them private and provide
236  * method that returns Region accessor. ?? Speed concerns? This introduce needless copies of mesh pointer since
237  * in most cases we do not need methods label, id, dim where whole RegionnDB is needed. Thus possibly make
238  * class RegionIdx (without Mesh) and make some asking methods for label, id in RegionDB.
239  * - in RegionDB make support for sets:
240 
241 /// map set name to lists of indices of its regions
242 typedef std::vector<RegionIdx> RegionSet;
243 std::map<std::string, RegionSet > sets_;
244 
245 /// Add region to given set. Creat the set if it does not exist.
246 add_to_set( const string& set_name, RegionIdx region);
247 /// Add a set into map, delete possible previous value, do not worry about slow copies of
248 /// the set array.
249 add_set( const string& set_name, const RegionSet & set);
250 RegionSet union( const string & set_name_1, const string & set_name_2); // sort + std::set_union
251 RegionSet intersection( const string & set_name_1, const string & set_name_2);
252 RegionSet difference( const string & set_name_1, const string & set_name_2);
253 const RegionSet & get_set(const string & set_name);
254 void read_sets_from_input(Input::Record rec); // see structure of RegionDB::region_set_input_type
255 
256 region_sets = [
257  { name="set name",
258  region_ids=[ int ...],
259  region_labels= ["..."], // these are merger together
260 
261  union=["set_1", "set_2"], // later overwrites previous
262  intersection=
263  difference=
264  }
265 ]
266  *
267  * - MODIFICATION IN MESH AND READER:
268  * Mesh reading proccess:
269  * 1) Read PhysicalNames form GMSH file, populate RegionDB (DONE in GMSH reader, may need small modifications)
270  * 2) Read region definitions from the input, see
271  *
272  * typedef std::map<unsigned int, unsigned int> MapElementIDToRegionID;
273  * RegionDB::read_regions_from_input(Input::Array region_list, MapElementIDToRegionID &map);
274  *
275  *
276  * (TODO in RegionDB, also creates (and return to Mesh) element regions modification map: std::map< unsigned int, RegionIdx>
277  * that maps element IDs to the new region names, GMSH reader should have setter method to accept this map
278  * and modify the elements during reading)
279  * 3) Read region sets - TODO in RegionDB
280  * 4) Read boundary key of the Mesh record and mark appropriate regions as boundary (TODO in RegionDB)
281  * 5) Read nodes (DONE in GMSH reader)
282  * 6) Read elements, per element:
283  * - possibly modify region according map
284  * - find element ID:
285  * if found: add_region(ID, get_label, dim, get_boundary_flag) // possibly set dimension of the region if it is undefined
286  * not found: add_region(ID, default_label, dim, false)
287  * - if region is boundary put element into Mesh::bc_elements
288  * else put it into Mesh::element
289  * ---
290  * 7) Setup topology - we has to connect Boundary with existing bc_elements, and add the remaining elements,
291  * after we remove support for old bCD files we may skip creation of remaining boundary elements since there will be no way how to set
292  * BC on them.
293  *
294  */
295 
296 class RegionDB {
297 public:
298  /**
299  * Map representing the relevance of elements to regions
300  */
302 
303  /**
304  * Format of input record which defined elements and their affiliation to region sets
305  */
307  /**
308  * Format of input record which defined regions and their affiliation to region sets
309  */
311 
312  TYPEDEF_ERR_INFO( EI_Label, const std::string);
313  TYPEDEF_ERR_INFO( EI_ID, unsigned int);
314  TYPEDEF_ERR_INFO( EI_IDOfOtherLabel, unsigned int);
315  TYPEDEF_ERR_INFO( EI_LabelOfOtherID, const std::string);
316  DECLARE_EXCEPTION( ExcAddingIntoClosed, << "Can not add label=" << EI_Label::qval << " into closed MaterialDispatch.\n");
317  //DECLARE_EXCEPTION( ExcSizeWhileOpen, << "Can not get size of MaterialDispatch yet open.");
318  DECLARE_EXCEPTION( ExcNonuniqueID, << "Non-unique ID during add of region id: " << EI_ID::val << ", label: " << EI_Label::qval << "\n" \
319  << "other region with same ID but different label: " << EI_LabelOfOtherID::qval << " already exists\n");
320  DECLARE_EXCEPTION( ExcNonuniqueLabel, << "Non-unique label during add of region id: " << EI_ID::val << ", label: " << EI_Label::qval << "\n" \
321  << "other region with same label but different ID: " << EI_IDOfOtherLabel::val << " already exists\n");
322  DECLARE_EXCEPTION( ExcInconsistentBoundary, << "Inconsistent add of region with id: " << EI_ID::val << ", label: " << EI_Label::qval << "\n" \
323  << "both ID and label match an existing region with different boundary flag.");
324  DECLARE_EXCEPTION( ExcInconsistentDimension, << "Inconsistent add of region with id: " << EI_ID::val << ", label: " << EI_Label::qval << "\n" \
325  << "both ID and label match an existing region with different dimension.");
326 
327  DECLARE_EXCEPTION( ExcCantAdd, << "Can not add new region into DB, id: " << EI_ID::val <<", label: " << EI_Label::qval);
328 
329  DECLARE_EXCEPTION( ExcUnknownSet, << "Operation with unknown region set: " << EI_Label::qval );
330 
331  DECLARE_INPUT_EXCEPTION( ExcUnknownSetOperand, << "Operation with unknown region set: " << EI_Label::qval);
332 
333  TYPEDEF_ERR_INFO( EI_NumOp, unsigned int);
334  DECLARE_INPUT_EXCEPTION( ExcWrongOpNumber, << "Wrong number of operands. Expect 2, given: " << EI_NumOp::val);
335 
336  /// Default constructor
337  RegionDB();
338 
339  /**
340  * Introduce an artificial limit to keep all material indexed arrays
341  * of reasonable size.
342  */
343  static const unsigned int max_n_regions = 64000;
344 
345  /// Undefined dimension for regions introduced from mesh input record.
346  /// Dimensions 0,1,2,3 are valid.
347  static const unsigned int undefined_dim = 10;
348 
349 
350  /**
351  * This method adds new region into the database and returns its index. This requires full
352  * specification of the region that is given in PhysicalNames section of the GMSH MSH format
353  * or in Mesh input record.
354  *
355  * If ID or label are found in the DB, we distinguish following cases:
356  * 1) ID is found, label is not found : warning ID has already assigned label
357  * 2) ID is not found, label is found : report error - assigning same label to different IDs
358  * 3) both ID and label are found, in same region : check remaining data, return existing region
359  * 4) , in different : warning ID has already assigned label
360  *
361  * Parameter @p id is any unique non-negative integer, parameter @p label is unique string identifier of the region,
362  * @p dim is dimension of reference elements in the region and @p boundary is true if the region consist of boundary elements
363  * (where one can apply boundary condition).
364  *
365  */
366  Region add_region(unsigned int id, const std::string &label, unsigned int dim, bool boundary);
367 
368  /**
369  * As the previous, but set the 'boundary; flag according to the label (labels starting with dot '.' are boundary).
370  * Used in read_regions_from_input ( with undefined dimension) to read regions given in 'regions' key of the 'mesh' input record.
371  */
372  Region add_region(unsigned int id, const std::string &label, unsigned int dim);
373 
374  /**
375  * As the previous, but generates automatic label of form 'region_ID' if the region with same ID is not already present. Set bulk region.
376  * Meant to be used when reading elements from MSH file. Again, if the region is defined already, we just check consistency.
377  */
378  Region add_region(unsigned int id, unsigned int dim);
379 
380  /**
381  * Returns a @p Region with given @p label. If it is not found it returns @p undefined Region.
382  */
383  Region find_label(const std::string &label) const;
384 
385  /**
386  * Returns a @p Region with given @p id. If it is not found it returns @p undefined Region.
387  */
388  Region find_id(unsigned int id) const;
389 
390  /**
391  * Return original label for given index @p idx.
392  */
393  const std::string &get_label(unsigned int idx) const;
394 
395  /**
396  * Return original ID for given index @p idx.
397  */
398  unsigned int get_id(unsigned int idx) const;
399 
400  /**
401  * Return dimension of region with given index @p idx.
402  */
403  unsigned int get_dim(unsigned int idx) const;
404 
405  /**
406  * Close this class for adding labels. This is necessary to return correct size
407  * for material indexed arrays and vectors. After calling this method you can
408  * call method @p size and method @p idx_of_label rise an exception for any unknown label.
409  */
410  void close();
411 
412  /**
413  * Returns maximal index + 1
414  */
415  unsigned int size() const;
416 
417  /**
418  * Returns total number boundary regions.
419  */
420  unsigned int boundary_size() const;
421 
422  /**
423  * Returns total number bulk regions.
424  */
425  unsigned int bulk_size() const;
426 
427  /**
428  * Returns implicit bulk region.
429  */
430  //Region implicit_bulk() const
431  //{ return implicit_bulk_; }
432 
433  /**
434  * Returns implicit boundary region. Is used for boundary elements created by Flow123d itself.
435  * This region has label "IMPLICIT_BOUNDARY".
436  */
438 
439  /*
440  * Add region to given set. Create the set if it does not exist.
441  *
442  * @param set_name Set to which it is added region
443  * @param region Added region
444  */
445  void add_to_set( const string& set_name, Region region);
446 
447  /**
448  * Add a set into map, delete possible previous value.
449  *
450  * @param set_name Name of added set
451  * @param set Added RegionSet
452  */
453  void add_set( const string& set_name, const RegionSet & set);
454 
455  /**
456  * Get RegionSets of specified names and create their union
457  *
458  * @param set_name_1 Name of first RegionSet
459  * @param set_name_2 Name of second RegionSet
460  * @return RegionSet created of union operation
461  */
462  RegionSet union_sets( const string & set_name_1, const string & set_name_2);
463 
464  /**
465  * Get RegionSets of specified names and create their intersection.
466  * Throws ExcUnknownSet for invalid name.
467  *
468  * @param set_name_1 Name of first RegionSet
469  * @param set_name_2 Name of second RegionSet
470  * @return RegionSet created of intersection operation
471  */
472  RegionSet intersection( const string & set_name_1, const string & set_name_2);
473 
474  /**
475  * Get RegionSets of specified names and create their difference
476  * Throws ExcUnknownSet for invalid name.
477  *
478  * @param set_name_1 Name of first RegionSet
479  * @param set_name_2 Name of second RegionSet
480  * @return RegionSet created of difference operation
481  */
482  RegionSet difference( const string & set_name_1, const string & set_name_2);
483 
484  /**
485  * Get region set of specified name. Three sets are defined by default:
486  * "ALL" - set of all regions both bulk and boundary.
487  * "BULK" - set of all bulk regions
488  * "BOUNDARY" - set of all boundary regions
489  *
490  * @param set_name Name of set
491  * @return RegionSet of specified name. Returns Empty vector if the set of given name doesn't exist.
492  */
493  RegionSet get_region_set(const string & set_name) const;
494 
495  /**
496  * Reads region sets defined by user in input file
497  * Format of input record is defined in variable RegionDB::region_set_input_type
498  *
499  * @param arr Array input records which define region sets
500  */
502 
503  /**
504  * Reads elements and their affiliation to region sets defined by user in input file
505  * Format of input record is defined in variable RegionDB::region_input_type
506  *
507  * @param region_list Array input records which define region sets and elements
508  * @param map Map to which is loaded data
509  */
511 
512 
513 private:
514 
515 
516 
517  /// One item in region database
518  struct RegionItem {
519  RegionItem(unsigned int index, unsigned int id, const std::string &label, unsigned int dim)
520  : index(index), id(id), label(label), dim_(dim) {}
521 
522  // unique identifiers
523  unsigned int index;
524  unsigned int id;
525  std::string label;
526  // data
527  unsigned int dim_;
528  };
529 
530  // tags
531  struct ID {};
532  struct Label {};
533  struct Index {};
534 
535  /// Region database
536  typedef BMI::multi_index_container<
537  RegionItem,
538  BMI::indexed_by<
539  // access by index
540  // Can not use random access without introducing "empty" RegionItems to fill holes in the case
541  // n_boundary != n_bulk. Empty items must be unsince we may have empty (and unmodifiable) holes ?? why)
542  //
543  // we need O(1) access
544  //BMI::random_access< BMI::tag<RandomIndex > >,
545  BMI::ordered_unique< BMI::tag<Index>, BMI::member<RegionItem, unsigned int, &RegionItem::index > >,
546  // use hashing for IDs, to get O(1) find complexity .. necessary for large meshes
547  BMI::hashed_unique< BMI::tag<ID>, BMI::member<RegionItem, unsigned int, &RegionItem::id> >,
548  // ordered access (like stl::map) by label
549  BMI::ordered_unique< BMI::tag<Label>, BMI::member<RegionItem, std::string, &RegionItem::label> >
550  >
552 
553  // ID and Label index iterators
554  typedef RegionTable::index<Label>::type::iterator LabelIter;
555  typedef RegionTable::index<ID>::type::iterator IDIter;
556 
557 
558  /// Consistency check in common use by various add_region methods.
559  void check_dim_consistency(IDIter it_id, unsigned int dim);
560 
561 
562  /// Database of all regions (both boundary and bulk).
564 
565  /// flag for closed database, no regions can be added, but you can add region sets
566  bool closed_;
567  /// Number of boundary regions
568  unsigned int n_boundary_;
569  /// Number of bulk regions
570  unsigned int n_bulk_;
571 
572  /// Map of region sets
574 
575  /// Make part of general RegionSet table.
577 
578  /**
579  * Implicit bulk and boundary regions. For GMSH mesh format only implicit_boundary is used for boundary elements
580  * that are not explicitly in the mesh file.
581  */
583 
584  /**
585  * Prepare region sets for union, intersection and difference operation.
586  * Get sets of names set_name_1 and set_name_2 and sort them.
587  * Throws ExcUnknownSet if the set with given name does not exist.
588  */
589  void prepare_sets( const string & set_name_1, const string & set_name_2, RegionSet & set_1, RegionSet & set_2);
590 
591  /**
592  * Read two operands from input array of strings and check if given names
593  * are existing sets. Return pair of checked set names.
594  */
595  pair<string,string> get_and_check_operands(const Input::Array & operands);
596 
597 };
598 
599 
600 
601 #endif /* REGION_HH_ */