Flow123d  master-eb60559f2
mh_dofhandler.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file mh_dofhandler.hh
15  * @brief
16  */
17 
18 #ifndef MH_DOFHANDLER_HH_
19 #define MH_DOFHANDLER_HH_
20 
21 #include <ext/alloc_traits.h> // for __alloc_traits<>::value_...
22 #include <sys/types.h> // for uint
23 #include <memory> // for shared_ptr, __shared_ptr
24 #include <unordered_map> // for unordered_map
25 #include <vector> // for vector
26 #include <armadillo>
27 #include "la/distribution.hh" // for Distribution
28 #include "system/index_types.hh" // for LongIdx
29 #include "mesh/accessors.hh" // for ElementAccessor, Side::edge_idx
30 #include "mesh/elements.h" // for Element::side, Element::dim
31 #include "mesh/mesh.h" // for Mesh
32 #include "mesh/region.hh" // for Region
33 #include "fem/dh_cell_accessor.hh" // for DHCellAccessor
34 
35 template <int spacedim> class LocalElementAccessorBase;
36 class DarcyMH;
37 class DarcyLMH;
38 class RichardsLMH;
39 
40 using namespace std;
41 
42 /// temporary solution to provide access to results
43 /// from DarcyFlowMH independent of mesh
45 public:
46  MH_DofHandler();
47  ~MH_DofHandler();
48  void reinit(Mesh *mesh);
49 
50  void prepare_parallel();
51  void make_row_numberings();
52 
53  void set_solution( double time, double * solution);
54 
55  inline double time_changed() const
56  { return time_; }
57 
58  unsigned int side_dof(const SideIter side) const;
59 
60  /// temporary replacement for DofHandler accessor, flux through given side
61  double side_flux(const Side &side) const;
62 
63  /// temporary replacement for DofHandler accessor, scalar (pressure) on edge of the side
64  double side_scalar(const Side &side) const;
65 
66  /// temporary replacement for DofHandler accessor, scalar (pressure) on element
67  double element_scalar( ElementAccessor<3> &ele ) const;
68 
69 protected:
71 
73  LongIdx *el_4_loc; //< array of idexes of local elements (in ordering matching the optimal global)
74  LongIdx *row_4_el; //< element index to matrix row
75  LongIdx *side_id_4_loc; //< array of ids of local sides
76  LongIdx *side_row_4_id; //< side id to matrix row
77  LongIdx *edge_4_loc; //< array of indexes of local edges
78  LongIdx *row_4_edge; //< edge index to matrix row
79 
80  // parallel
81  Distribution *edge_ds; //< optimal distribution of edges
82  Distribution *el_ds; //< optimal distribution of elements
83  Distribution *side_ds; //< optimal distribution of elements
84 
85 
86  /// Maps mesh index of the edge to the edge index in the mesh portion local to the processor.
87  /// Temporary solution until we have parallel mesh which should provide such information.
88  std::unordered_map<unsigned int, unsigned int> edge_new_local_4_mesh_idx_;
89 
90  double * mh_solution;
91  double time_;
92 
94  friend DarcyMH;
95  friend DarcyLMH;
96  friend RichardsLMH;
97 };
98 
99 
100 
101 typedef unsigned int uint;
102 
103 template <int spacedim>
105 public:
106 
108  : dh_cell_(dh_cell), global_indices_(dh_cell_.dh()->max_elem_dofs())
109  {
110  n_indices_ = dh_cell_.get_dof_indices(global_indices_);
111  }
112 
113  inline DHCellAccessor dh_cell() const {
114  return dh_cell_;
115  }
116 
117  uint dim() const {
118  return dh_cell_.dim();
119  }
120 
121  uint n_sides() const {
122  return element_accessor()->n_sides();
123  }
124 
126  return dh_cell_.elm();
127  }
128 
129  const arma::vec3 centre() const {
130  return element_accessor().centre();
131  }
132 
133  double measure() const {
134  return element_accessor().measure();
135  }
136 
137  Region region() const {
138  return element_accessor().region();
139  }
140 
142  return element_accessor().idx();
143  }
144 
145  uint ele_local_idx() const {
146  return dh_cell_.local_idx();
147  }
148 
150  return global_indices_[n_indices_/2];
151  }
152 
154  return dh_cell_.get_loc_dof_indices()[n_indices_/2];
155  }
156 
158  return global_indices_[(n_indices_+1)/2+i];
159  }
160 
162  return dh_cell_.get_loc_dof_indices()[(n_indices_+1)/2+i];
163  }
164 
166  return element_accessor().side(i);
167  }
168 
170  return global_indices_[i];
171  }
172 
174  return dh_cell_.get_loc_dof_indices()[i];
175  }
176 
177 private:
181 };
182 
183 /**
184  * This is prototype of further much more complex and general accessor templated by
185  * element dimension. In fact we shall need an accessor for every kind of element interaction integral.
186  */
187 template <int spacedim, int dim>
189 public:
191  : LocalElementAccessorBase<spacedim>(dh, loc_ele_idx)
192  {}
193 };
194 
195 
196 
197 #endif /* MH_DOFHANDLER_HH_ */
LocalElementAccessorBase::edge_local_row
uint edge_local_row(uint i)
Definition: mh_dofhandler.hh:161
MH_DofHandler::edge_4_loc
LongIdx * edge_4_loc
Definition: mh_dofhandler.hh:77
LocalElementAccessorBase::dh_cell_
DHCellAccessor dh_cell_
Definition: mh_dofhandler.hh:178
LocalElementAccessorBase::ele_global_idx
uint ele_global_idx()
Definition: mh_dofhandler.hh:141
distribution.hh
Support classes for parallel programing.
MH_DofHandler::RichardsLMH
friend RichardsLMH
Definition: mh_dofhandler.hh:96
MH_DofHandler::edge_new_local_4_mesh_idx_
std::unordered_map< unsigned int, unsigned int > edge_new_local_4_mesh_idx_
Definition: mh_dofhandler.hh:88
MH_DofHandler::mesh_
Mesh * mesh_
Definition: mh_dofhandler.hh:72
std::vector
Definition: doxy_dummy_defs.hh:7
ElementAccessor< 3 >
LocalElementAccessorBase
Definition: mh_dofhandler.hh:35
arma::vec3
Definition: doxy_dummy_defs.hh:17
LocalElementAccessorBase::LocalElementAccessorBase
LocalElementAccessorBase(DHCellAccessor dh_cell)
Definition: mh_dofhandler.hh:107
MH_DofHandler::row_4_el
LongIdx * row_4_el
Definition: mh_dofhandler.hh:74
uint
unsigned int uint
Definition: mh_dofhandler.hh:101
LocalElementAccessorBase::measure
double measure() const
Definition: mh_dofhandler.hh:133
index_types.hh
LocalElementAccessorBase::dh_cell
DHCellAccessor dh_cell() const
Definition: mh_dofhandler.hh:113
MH_DofHandler::edge_ds
Distribution * edge_ds
Definition: mh_dofhandler.hh:81
MH_DofHandler::DarcyLMH
friend DarcyLMH
Definition: mh_dofhandler.hh:95
LocalElementAccessorBase::centre
const arma::vec3 centre() const
Definition: mh_dofhandler.hh:129
Region
Definition: region.hh:145
MH_DofHandler::elem_side_to_global
vector< vector< unsigned int > > elem_side_to_global
Definition: mh_dofhandler.hh:70
LocalElementAccessor
Definition: mh_dofhandler.hh:188
Distribution
Definition: distribution.hh:50
dh_cell_accessor.hh
accessors.hh
LocalElementAccessorBase::ele_row
uint ele_row()
Definition: mh_dofhandler.hh:149
elements.h
Side
Definition: accessors.hh:390
MH_DofHandler::side_ds
Distribution * side_ds
Definition: mh_dofhandler.hh:83
LocalElementAccessorBase::global_indices_
std::vector< LongIdx > global_indices_
Definition: mh_dofhandler.hh:179
mesh.h
MH_DofHandler::side_id_4_loc
LongIdx * side_id_4_loc
Definition: mh_dofhandler.hh:75
LocalElementAccessor::LocalElementAccessor
LocalElementAccessor(MH_DofHandler &dh, uint loc_ele_idx)
Definition: mh_dofhandler.hh:190
MH_DofHandler::el_ds
Distribution * el_ds
Definition: mh_dofhandler.hh:82
LongIdx
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
LocalElementAccessorBase::side_row
uint side_row(uint i)
Definition: mh_dofhandler.hh:169
LocalElementAccessorBase::edge_row
uint edge_row(uint i)
Definition: mh_dofhandler.hh:157
Mesh
Definition: mesh.h:361
MH_DofHandler::mh_solution
double * mh_solution
Definition: mh_dofhandler.hh:90
DHCellAccessor
Cell accessor allow iterate over DOF handler cells.
Definition: dh_cell_accessor.hh:43
LocalElementAccessorBase::dim
uint dim() const
Definition: mh_dofhandler.hh:117
std
Definition: doxy_dummy_defs.hh:5
DarcyLMH
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
Definition: darcy_flow_lmh.hh:133
MH_DofHandler::time_
double time_
Definition: mh_dofhandler.hh:91
LocalElementAccessorBase::side_local_row
uint side_local_row(uint i)
Definition: mh_dofhandler.hh:173
region.hh
RichardsLMH
Edge lumped mixed-hybrid solution of unsteady Darcy flow.
Definition: richards_lmh.hh:65
LocalElementAccessorBase::element_accessor
ElementAccessor< 3 > element_accessor() const
Definition: mh_dofhandler.hh:125
LocalElementAccessorBase::ele_local_row
uint ele_local_row()
Definition: mh_dofhandler.hh:153
MH_DofHandler::DarcyMH
friend DarcyMH
Definition: mh_dofhandler.hh:94
MH_DofHandler::time_changed
double time_changed() const
Definition: mh_dofhandler.hh:55
LocalElementAccessorBase::n_sides
uint n_sides() const
Definition: mh_dofhandler.hh:121
LocalElementAccessorBase::region
Region region() const
Definition: mh_dofhandler.hh:137
LocalElementAccessorBase::ele_local_idx
uint ele_local_idx() const
Definition: mh_dofhandler.hh:145
DarcyMH
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
Definition: darcy_flow_mh.hh:125
MH_DofHandler::el_4_loc
LongIdx * el_4_loc
Definition: mh_dofhandler.hh:73
SideIter
Definition: accessors.hh:490
MH_DofHandler
Definition: mh_dofhandler.hh:44
MH_DofHandler::row_4_edge
LongIdx * row_4_edge
Definition: mh_dofhandler.hh:78
LocalElementAccessorBase::n_indices_
uint n_indices_
Definition: mh_dofhandler.hh:180
LocalElementAccessorBase::side
SideIter side(uint i)
Definition: mh_dofhandler.hh:165
MH_DofHandler::side_row_4_id
LongIdx * side_row_4_id
Definition: mh_dofhandler.hh:76