Flow123d  JS_before_hm-929-gaeebe69
assembly_mh.hh
Go to the documentation of this file.
1 /*
2  * darcy_flow_assembly.hh
3  *
4  * Created on: Apr 21, 2016
5  * Author: jb
6  */
7 
8 #ifndef SRC_FLOW_DARCY_FLOW_ASSEMBLY_HH_
9 #define SRC_FLOW_DARCY_FLOW_ASSEMBLY_HH_
10 
11 #include "system/index_types.hh"
12 #include "mesh/mesh.h"
13 #include "mesh/accessors.hh"
14 #include "mesh/neighbours.h"
15 #include "fem/fe_p.hh"
16 #include "fem/fe_values.hh"
17 #include "fem/fe_rt.hh"
18 #include "fem/fe_values_views.hh"
20 
21 #include "la/linsys.hh"
22 #include "la/linsys_PETSC.hh"
23 #include "la/linsys_BDDC.hh"
24 #include "la/schur.hh"
25 
26 #include "la/local_system.hh"
27 
28 #include "coupling/balance.hh"
29 #include "flow/darcy_flow_mh.hh"
30 #include "flow/mortar_assembly.hh"
31 
32 
33 /** Common abstract class for the assembly routines in Darcy flow.
34  * Is implemented in DarcyMH, DarcyLMH and RichardsLMH assembly classes,
35  * which are independent of each other.
36  */
38 {
39 public:
40  virtual void fix_velocity(const DHCellAccessor& dh_cell) = 0;
41  virtual void assemble(const DHCellAccessor& dh_cell) = 0;
42  virtual void assemble_reconstruct(const DHCellAccessor& dh_cell) = 0;
43 
44  /// Updates water content in Richards.
45  virtual void update_water_content(const DHCellAccessor& dh_cell) = 0;
46 
47  /**
48  * Generic creator of multidimensional assembly, i.e. vector of
49  * particular assembly objects.
50  */
51  template< template<int dim> class Impl, class Data>
52  static MultidimAssembly create(Data data) {
53  return { std::make_shared<Impl<1> >(data),
54  std::make_shared<Impl<2> >(data),
55  std::make_shared<Impl<3> >(data) };
56  }
57 
58  virtual ~AssemblyBase() {}
59 };
60 
61 
62 template <int dim>
64 private:
65  // assembly face integrals (BC)
68 public:
70  : side_quad_(dim, 1),
71  fe_p_disc_(0)
72  {
73  fe_side_values_.initialize(side_quad_, fe_p_disc_, update_normal_vectors);
74  }
76 
77 };
78 
79 
80 /** MH version of Darcy flow assembly. It is supposed not to be improved anymore,
81  * however it is kept functioning aside of the LMH lumped version until
82  * the LMH version is stable and optimized.
83  */
84 template<int dim>
85 class AssemblyMH : public AssemblyBase
86 {
87 public:
88  typedef std::shared_ptr<DarcyMH::EqData> AssemblyDataPtrMH;
89 
90  AssemblyMH<dim>(AssemblyDataPtrMH data)
91  : quad_(dim, 3),
92  velocity_interpolation_quad_(dim, 0), // veloctiy values in barycenter
93  ad_(data),
94  loc_system_(size(), size()),
95  loc_system_vb_(2,2)
96 
97  {
98  fe_values_.initialize(quad_, fe_rt_,
100  velocity_interpolation_fv_.initialize(velocity_interpolation_quad_, fe_rt_, update_values | update_quadrature_points);
101 
102  // local numbering of dofs for MH system
103  unsigned int nsides = dim+1;
104  loc_side_dofs.resize(nsides);
105  loc_ele_dof = nsides;
106  loc_edge_dofs.resize(nsides);
107  for(unsigned int i = 0; i < nsides; i++){
108  loc_side_dofs[i] = i;
109  loc_edge_dofs[i] = nsides + i + 1;
110  }
111  //DebugOut() << print_var(this) << print_var(side_quad_.size());
112 
113  // create local sparsity pattern
114  arma::umat sp(size(),size());
115  sp.zeros();
116  sp.submat(0, 0, nsides, nsides).ones();
117  sp.diag().ones();
118  // armadillo 8.4.3 bug with negative sub diagonal index
119  // sp.diag(nsides+1).ones();
120  // sp.diag(-nsides-1).ones();
121  // sp.print();
122 
123  sp.submat(0, nsides+1, nsides-1, size()-1).diag().ones();
124  sp.submat(nsides+1, 0, size()-1, nsides-1).diag().ones();
125 
126  loc_system_.set_sparsity(sp);
127 
128  // local system 2x2 for vb neighbourings is full matrix
129  // this matrix cannot be influenced by any BC (no elimination can take place)
130  sp.ones(2,2);
131  loc_system_vb_.set_sparsity(sp);
132 
133  if (ad_->mortar_method_ == DarcyFlowInterface::MortarP0) {
134  mortar_assembly = std::make_shared<P0_CouplingAssembler>(ad_);
135  } else if (ad_->mortar_method_ == DarcyFlowInterface::MortarP1) {
136  mortar_assembly = std::make_shared<P1_CouplingAssembler>(ad_);
137  }
138 
139  }
140 
141  void assemble_reconstruct(const DHCellAccessor&) override
142  {};
143  void update_water_content(const DHCellAccessor&) override
144  {};
145 
146  ~AssemblyMH<dim>() override
147  {}
148 
149 // LocalSystem& get_local_system() override
150 // { return loc_system_;}
151 
152  void fix_velocity(const DHCellAccessor& dh_cell) override
153  {
154  if (mortar_assembly)
155  mortar_assembly->fix_velocity(dh_cell);
156  }
157 
158  void assemble(const DHCellAccessor& dh_cell) override
159  {
160  ASSERT_EQ_DBG(dh_cell.dim(), dim);
161  loc_system_.reset();
162 
163  set_dofs_and_bc(dh_cell);
164 
165  assemble_sides(dh_cell);
166  assemble_element(dh_cell);
167 
168  loc_system_.eliminate_solution();
169  ad_->lin_sys->set_local_system(loc_system_);
170 
171  assembly_dim_connections(dh_cell);
172 
173  if (ad_->balance != nullptr)
174  add_fluxes_in_balance_matrix(dh_cell);
175 
176  if (mortar_assembly)
177  mortar_assembly->assembly(dh_cell);
178  }
179 
180  void assembly_local_vb(ElementAccessor<3> ele, DHCellSide neighb_side) //override
181  {
182  ASSERT_LT_DBG(ele->dim(), 3);
183  //DebugOut() << "alv " << print_var(this);
184  //START_TIMER("Assembly<dim>::assembly_local_vb");
185  // compute normal vector to side
186  arma::vec3 nv;
187  ElementAccessor<3> ele_higher = ad_->mesh->element_accessor( neighb_side.element().idx() );
188  ngh_values_.fe_side_values_.reinit(neighb_side.side());
189  nv = ngh_values_.fe_side_values_.normal_vector(0);
190 
191  double value = ad_->sigma.value( ele.centre(), ele) *
192  2*ad_->conductivity.value( ele.centre(), ele) *
193  arma::dot(ad_->anisotropy.value( ele.centre(), ele)*nv, nv) *
194  ad_->cross_section.value( neighb_side.centre(), ele_higher ) * // cross-section of higher dim. (2d)
195  ad_->cross_section.value( neighb_side.centre(), ele_higher ) /
196  ad_->cross_section.value( ele.centre(), ele ) * // crossection of lower dim.
197  neighb_side.measure();
198 
199  loc_system_vb_.add_value(0,0, -value);
200  loc_system_vb_.add_value(0,1, value);
201  loc_system_vb_.add_value(1,0, value);
202  loc_system_vb_.add_value(1,1, -value);
203  }
204 
205 protected:
206  static unsigned int size()
207  {
208  // dofs: velocity, pressure, edge pressure
210  }
211 
212  void set_dofs_and_bc(const DHCellAccessor& dh_cell){
213 
214  local_dofs_ = dh_cell.get_loc_dof_indices();
215  global_dofs_.resize(dh_cell.n_dofs());
216  dh_cell.get_dof_indices(global_dofs_);
217 
218  const ElementAccessor<3> ele = dh_cell.elm();
219 
220  //set global dof for element (pressure)
221  loc_system_.row_dofs[loc_ele_dof] = loc_system_.col_dofs[loc_ele_dof] = global_dofs_[loc_ele_dof];
222 
223  //shortcuts
224  const unsigned int nsides = ele->n_sides();
225  LinSys *ls = ad_->lin_sys;
226 
227  unsigned int side_row, edge_row;
228 
229  dirichlet_edge.resize(nsides);
230  for (unsigned int i = 0; i < nsides; i++) {
231 
232  side_row = loc_side_dofs[i]; //local
233  edge_row = loc_edge_dofs[i]; //local
234  loc_system_.row_dofs[side_row] = loc_system_.col_dofs[side_row] = global_dofs_[side_row]; //global
235  loc_system_.row_dofs[edge_row] = loc_system_.col_dofs[edge_row] = global_dofs_[edge_row]; //global
236 
237  dirichlet_edge[i] = 0;
238  Side side = *dh_cell.elm().side(i);
239  if (side.is_boundary()) {
240  Boundary bcd = side.cond();
241  ElementAccessor<3> b_ele = bcd.element_accessor();
242  DarcyMH::EqData::BC_Type type = (DarcyMH::EqData::BC_Type)ad_->bc_type.value(b_ele.centre(), b_ele);
243 
244  double cross_section = ad_->cross_section.value(ele.centre(), ele);
245 
246  if ( type == DarcyMH::EqData::none) {
247  // homogeneous neumann
248  } else if ( type == DarcyMH::EqData::dirichlet ) {
249  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
250  loc_system_.set_solution(loc_edge_dofs[i],bc_pressure,-1);
251  dirichlet_edge[i] = 1;
252 
253  } else if ( type == DarcyMH::EqData::total_flux) {
254  // internally we work with outward flux
255  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
256  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
257  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
258 
259  dirichlet_edge[i] = 2; // to be skipped in LMH source assembly
260  loc_system_.add_value(edge_row, edge_row,
261  -b_ele.measure() * bc_sigma * cross_section,
262  (bc_flux - bc_sigma * bc_pressure) * b_ele.measure() * cross_section);
263  }
264  else if (type==DarcyMH::EqData::seepage) {
265  ad_->is_linear=false;
266 
267  unsigned int loc_edge_idx = bcd.bc_ele_idx();
268  char & switch_dirichlet = ad_->bc_switch_dirichlet[loc_edge_idx];
269  double bc_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
270  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
271  double side_flux = bc_flux * b_ele.measure() * cross_section;
272 
273  // ** Update BC type. **
274  if (switch_dirichlet) {
275  // check and possibly switch to flux BC
276  // The switch raise error on the corresponding edge row.
277  // Magnitude of the error is abs(solution_flux - side_flux).
278  ASSERT_DBG(ad_->dh_->distr()->is_local(global_dofs_[side_row]))(global_dofs_[side_row]);
279  unsigned int loc_side_row = local_dofs_[side_row];
280  double & solution_flux = ls->get_solution_array()[loc_side_row];
281 
282  if ( solution_flux < side_flux) {
283  //DebugOut().fmt("x: {}, to neum, p: {} f: {} -> f: {}\n", b_ele.centre()[0], bc_pressure, solution_flux, side_flux);
284  solution_flux = side_flux;
285  switch_dirichlet=0;
286  }
287  } else {
288  // check and possibly switch to pressure BC
289  // TODO: What is the appropriate DOF in not local?
290  // The switch raise error on the corresponding side row.
291  // Magnitude of the error is abs(solution_head - bc_pressure)
292  // Since usually K is very large, this error would be much
293  // higher then error caused by the inverse switch, this
294  // cause that a solution with the flux violating the
295  // flux inequality leading may be accepted, while the error
296  // in pressure inequality is always satisfied.
297  ASSERT_DBG(ad_->dh_->distr()->is_local(global_dofs_[edge_row]))(global_dofs_[edge_row]);
298  unsigned int loc_edge_row = local_dofs_[edge_row];
299  double & solution_head = ls->get_solution_array()[loc_edge_row];
300 
301  if ( solution_head > bc_pressure) {
302  //DebugOut().fmt("x: {}, to dirich, p: {} -> p: {} f: {}\n",b_ele.centre()[0], solution_head, bc_pressure, bc_flux);
303  solution_head = bc_pressure;
304  switch_dirichlet=1;
305  }
306  }
307 
308  // ** Apply BCUpdate BC type. **
309  // Force Dirichlet type during the first iteration of the unsteady case.
310  if (switch_dirichlet || ad_->force_no_neumann_bc ) {
311  //DebugOut().fmt("x: {}, dirich, bcp: {}\n", b_ele.centre()[0], bc_pressure);
312  loc_system_.set_solution(loc_edge_dofs[i],bc_pressure, -1);
313  dirichlet_edge[i] = 1;
314  } else {
315  //DebugOut()("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], side_flux, bc_flux);
316  loc_system_.add_value(edge_row, side_flux);
317  }
318 
319  } else if (type==DarcyMH::EqData::river) {
320  ad_->is_linear=false;
321 
322  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
323  double bc_switch_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
324  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
325  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
326  ASSERT_DBG(ad_->dh_->distr()->is_local(global_dofs_[edge_row]))(global_dofs_[edge_row]);
327  unsigned int loc_edge_row = local_dofs_[edge_row];
328  double & solution_head = ls->get_solution_array()[loc_edge_row];
329 
330  // Force Robin type during the first iteration of the unsteady case.
331  if (solution_head > bc_switch_pressure || ad_->force_no_neumann_bc) {
332  // Robin BC
333  //DebugOut().fmt("x: {}, robin, bcp: {}\n", b_ele.centre()[0], bc_pressure);
334  loc_system_.add_value(edge_row, edge_row,
335  -b_ele.measure() * bc_sigma * cross_section,
336  b_ele.measure() * cross_section * (bc_flux - bc_sigma * bc_pressure) );
337  } else {
338  // Neumann BC
339  //DebugOut().fmt("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], bc_switch_pressure, bc_pressure);
340  double bc_total_flux = bc_flux + bc_sigma*(bc_switch_pressure - bc_pressure);
341 
342  loc_system_.add_value(edge_row, bc_total_flux * b_ele.measure() * cross_section);
343  }
344  }
345  else {
346  xprintf(UsrErr, "BC type not supported.\n");
347  }
348  }
349  loc_system_.add_value(side_row, edge_row, 1.0);
350  loc_system_.add_value(edge_row, side_row, 1.0);
351  }
352  }
353 
354  void assemble_sides(const DHCellAccessor& dh_cell)
355  {
356  const ElementAccessor<3> ele = dh_cell.elm();
357  double cs = ad_->cross_section.value(ele.centre(), ele);
358  double conduct = ad_->conductivity.value(ele.centre(), ele);
359  double scale = 1 / cs /conduct;
360 
361  assemble_sides_scale(dh_cell, scale);
362  }
363 
364  void assemble_sides_scale(const DHCellAccessor& dh_cell, double scale)
365  {
366  arma::vec3 &gravity_vec = ad_->gravity_vec_;
367  const ElementAccessor<3> ele = dh_cell.elm();
368 
369  fe_values_.reinit(ele);
370  unsigned int ndofs = fe_values_.n_dofs();
371  unsigned int qsize = fe_values_.n_points();
372  auto velocity = fe_values_.vector_view(0);
373 
374  for (unsigned int k=0; k<qsize; k++)
375  for (unsigned int i=0; i<ndofs; i++){
376  double rhs_val =
377  arma::dot(gravity_vec,velocity.value(i,k))
378  * fe_values_.JxW(k);
379  loc_system_.add_value(i, rhs_val);
380 
381  for (unsigned int j=0; j<ndofs; j++){
382  double mat_val =
383  arma::dot(velocity.value(i,k), //TODO: compute anisotropy before
384  (ad_->anisotropy.value(ele.centre(), ele )).i()
385  * velocity.value(j,k))
386  * scale * fe_values_.JxW(k);
387 
388  loc_system_.add_value(i, j, mat_val);
389  }
390  }
391 
392  // assemble matrix for weights in BDDCML
393  // approximation to diagonal of
394  // S = -C - B*inv(A)*B'
395  // as
396  // diag(S) ~ - diag(C) - 1./diag(A)
397  // the weights form a partition of unity to average a discontinuous solution from neighbouring subdomains
398  // to a continuous one
399  // it is important to scale the effect - if conductivity is low for one subdomain and high for the other,
400  // trust more the one with low conductivity - it will be closer to the truth than an arithmetic average
401  if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
402  const arma::mat& local_matrix = loc_system_.get_matrix();
403  for(unsigned int i=0; i < ndofs; i++) {
404  double val_side = local_matrix(i,i);
405  double val_edge = -1./local_matrix(i,i);
406 
407  unsigned int side_row = loc_system_.row_dofs[loc_side_dofs[i]];
408  unsigned int edge_row = loc_system_.row_dofs[loc_edge_dofs[i]];
409  static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( side_row, val_side );
410  static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( edge_row, val_edge );
411  }
412  }
413  }
414 
415 
417  // set block B, B': element-side, side-element
418 
419  for(unsigned int side = 0; side < loc_side_dofs.size(); side++){
420  loc_system_.add_value(loc_ele_dof, loc_side_dofs[side], -1.0);
421  loc_system_.add_value(loc_side_dofs[side], loc_ele_dof, -1.0);
422  }
423 
424  if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
425  double val_ele = 1.;
426  static_cast<LinSys_BDDC*>(ad_->lin_sys)->
427  diagonal_weights_set_value( loc_system_.row_dofs[loc_ele_dof], val_ele );
428  }
429  }
430 
431 
433  //D, E',E block: compatible connections: element-edge
434  auto ele = dh_cell.elm(); //ElementAccessor<3>
435 
436  // no Neighbours => nothing to asssemble here
437  if(dh_cell.elm()->n_neighs_vb() == 0) return;
438 
439  std::vector<LongIdx> higher_dim_dofs;
440  //DebugOut() << "adc " << print_var(this) << print_var(side_quad_.size());
441  unsigned int i = 0;
442  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
443  // every compatible connection adds a 2x2 matrix involving
444  // current element pressure and a connected edge pressure
445 
446  loc_system_vb_.reset();
447  loc_system_vb_.row_dofs[0] = loc_system_vb_.col_dofs[0] = global_dofs_[loc_ele_dof];
448 
449  // // read neighbor dofs (dh dofhandler)
450  // // neighbor cell owning neighb_side
451  DHCellAccessor dh_neighb_cell = neighb_side.cell();
452 
453  higher_dim_dofs.resize(dh_neighb_cell.n_dofs());
454  dh_neighb_cell.get_dof_indices(higher_dim_dofs);
455 
456  // local index of pedge dof on neighboring cell
457  // (dim+2) is number of edges of higher dim element
458  // TODO: replace with DHCell getter when available for FESystem component
459  const unsigned int t = dh_neighb_cell.n_dofs() - (dim+2) + neighb_side.side().side_idx();
460  loc_system_vb_.row_dofs[1] = loc_system_vb_.col_dofs[1] = higher_dim_dofs[t];
461 
462  assembly_local_vb(ele, neighb_side);
463 
464  loc_system_vb_.eliminate_solution();
465  ad_->lin_sys->set_local_system(loc_system_vb_);
466 
467  // update matrix for weights in BDDCML
468  if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
469  int ind = loc_system_vb_.row_dofs[1];
470  // there is -value on diagonal in block C!
471  double new_val = loc_system_vb_.get_matrix()(0,0);
472  static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( ind, new_val );
473  }
474  ++i;
475  }
476  }
477 
479 
480  for(DHCellSide dh_side : dh_cell.side_range()){
481  unsigned int sidx = dh_side.side_idx();
482  if (dh_side.side().is_boundary()) {
483  ad_->balance->add_flux_values(ad_->water_balance_idx, dh_side,
484  {local_dofs_[loc_side_dofs[sidx]]},
485  {1}, 0);
486  }
487  }
488  }
489 
490 
491 
492  // assembly volume integrals
496 
498 
499  // Interpolation of velocity into barycenters
502 
503  // data shared by assemblers of different dimension
504  AssemblyDataPtrMH ad_;
506 
511  unsigned int loc_ele_dof;
512 
513  std::shared_ptr<MortarAssemblyBase> mortar_assembly;
514 
517 };
518 
519 
520 #endif /* SRC_FLOW_DARCY_FLOW_ASSEMBLY_HH_ */
Shape function values.
Definition: update_flags.hh:87
Range< DHCellSide > side_range() const
Returns range of cell sides.
unsigned int loc_ele_dof
Definition: assembly_mh.hh:511
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
virtual void assemble(const DHCellAccessor &dh_cell)=0
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
std::vector< LongIdx > global_dofs_
Definition: assembly_mh.hh:515
QGauss velocity_interpolation_quad_
Definition: assembly_mh.hh:500
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
virtual void fix_velocity(const DHCellAccessor &dh_cell)=0
void assembly_dim_connections(const DHCellAccessor &dh_cell)
Definition: assembly_mh.hh:432
void assemble_sides_scale(const DHCellAccessor &dh_cell, double scale)
Definition: assembly_mh.hh:364
void assemble(const DHCellAccessor &dh_cell) override
Definition: assembly_mh.hh:158
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
Wrappers for linear systems based on MPIAIJ and MATIS format.
static unsigned int size()
Definition: assembly_mh.hh:206
unsigned int dim() const
Return dimension of element appropriate to cell.
Cell accessor allow iterate over DOF handler cells.
void set_dofs_and_bc(const DHCellAccessor &dh_cell)
Definition: assembly_mh.hh:212
Class FEValues calculates finite element data on the actual cells such as shape function values...
uint bc_ele_idx()
Definition: accessors.hh:345
LocalSystem loc_system_vb_
Definition: assembly_mh.hh:508
virtual void assemble_reconstruct(const DHCellAccessor &dh_cell)=0
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
double measure() const
void resize(unsigned int n_q_points)
Modify the number of quadrature points.
Definition: quadrature.hh:80
SideIter side(const unsigned int loc_index)
std::vector< unsigned int > loc_edge_dofs
Definition: assembly_mh.hh:510
std::vector< unsigned int > dirichlet_edge
Definition: assembly_mh.hh:505
unsigned int n_dofs() const
Return number of dofs on given cell.
Side side() const
Return Side of given cell and side_idx.
FEValues< 3 > velocity_interpolation_fv_
Definition: assembly_mh.hh:501
static constexpr bool value
Definition: json.hpp:87
std::shared_ptr< DarcyMH::EqData > AssemblyDataPtrMH
Definition: assembly_mh.hh:88
Symmetric Gauss-Legendre quadrature formulae on simplices.
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
unsigned int dim() const
Definition: elements.h:121
Transformed quadrature points.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
void assemble_sides(const DHCellAccessor &dh_cell)
Definition: assembly_mh.hh:354
ArmaMat< double, N, M > mat
Definition: armor.hh:864
virtual void update_water_content(const DHCellAccessor &dh_cell)=0
Updates water content in Richards.
bool is_boundary() const
Returns true for side on the boundary.
double * get_solution_array()
Definition: linsys.hh:325
FEValues< 3 > fe_values_
Definition: assembly_mh.hh:495
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:60
Definitions of basic Lagrangean finite elements with polynomial shape functions.
ElementAccessor< 3 > element() const
unsigned int n_sides() const
Definition: elements.h:132
#define xprintf(...)
Definition: system.hh:93
RangeConvert< DHNeighbSide, DHCellSide > neighb_sides() const
Returns range of neighbour cell of lower dimension corresponding to cell of higher dimension...
FEValues< 3 > fe_side_values_
Definition: assembly_mh.hh:75
LocalSystem loc_system_
Definition: assembly_mh.hh:507
double measure() const
Computes the measure of the element.
Shape function gradients.
Definition: update_flags.hh:95
Normal vectors.
Definition: system.hh:65
virtual ~AssemblyBase()
Definition: assembly_mh.hh:58
LocDofVec local_dofs_
Definition: assembly_mh.hh:516
void update_water_content(const DHCellAccessor &) override
Updates water content in Richards.
Definition: assembly_mh.hh:143
QGauss quad_
Definition: assembly_mh.hh:494
std::shared_ptr< MortarAssemblyBase > mortar_assembly
Definition: assembly_mh.hh:513
AssemblyDataPtrMH ad_
Definition: assembly_mh.hh:504
#define ASSERT_DBG(expr)
NeighSideValues< (dim< 3)?dim:2 > ngh_values_
Definition: assembly_mh.hh:497
Definitions of particular quadrature rules on simplices.
void assemble_reconstruct(const DHCellAccessor &) override
Definition: assembly_mh.hh:141
std::vector< unsigned int > loc_side_dofs
Definition: assembly_mh.hh:509
BasicData Data
Definition: format.h:848
void add_fluxes_in_balance_matrix(const DHCellAccessor &dh_cell)
Definition: assembly_mh.hh:478
mixed-hybrid model of linear Darcy flow, possibly unsteady.
static MultidimAssembly create(Data data)
Definition: assembly_mh.hh:52
void assemble_element(const DHCellAccessor &)
Definition: assembly_mh.hh:416
void assembly_local_vb(ElementAccessor< 3 > ele, DHCellSide neighb_side)
Definition: assembly_mh.hh:180
Boundary cond() const
FE_P_disc< dim+1 > fe_p_disc_
Definition: assembly_mh.hh:67
unsigned int n_neighs_vb() const
Return number of neighbours.
Definition: elements.h:67
FE_RT0< dim > fe_rt_
Definition: assembly_mh.hh:493
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
void fix_velocity(const DHCellAccessor &dh_cell) override
Definition: assembly_mh.hh:152
Quadrature * quad_
Quadrature used in assembling methods.
Definition: assembly_dg.hh:394
Definitions of Raviart-Thomas finite elements.
Side accessor allows to iterate over sides of DOF handler cell.
ElementAccessor< 3 > element_accessor()
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
Transformed quadrature weights.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
arma::vec3 centre() const
Side centre.