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