Flow123d  release_2.2.0-914-gf1a3a4f
darcy_flow_assembly.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 "mesh/mesh.h"
12 #include "fem/mapping_p1.hh"
13 #include "fem/fe_p.hh"
14 #include "fem/fe_values.hh"
15 #include "fem/fe_rt.hh"
16 #include "fem/fe_values_views.hh"
18 #include "flow/mh_dofhandler.hh"
19 
20 #include "la/linsys.hh"
21 #include "la/linsys_PETSC.hh"
22 #include "la/linsys_BDDC.hh"
23 #include "la/schur.hh"
24 
25 #include "la/local_system.hh"
26 
27 #include "coupling/balance.hh"
28 #include "flow/mortar_assembly.hh"
29 
30 
32 {
33 public:
34  typedef std::shared_ptr<DarcyMH::EqData> AssemblyDataPtr;
36 
37  virtual ~AssemblyBase() {}
38 
39  /**
40  * Generic creator of multidimensional assembly, i.e. vector of
41  * particular assembly objects.
42  */
43  template< template<int dim> class Impl >
44  static MultidimAssembly create(typename Impl<1>::AssemblyDataPtr data) {
45  return { std::make_shared<Impl<1> >(data),
46  std::make_shared<Impl<2> >(data),
47  std::make_shared<Impl<3> >(data) };
48 
49  }
50 
51 // virtual LocalSystem & get_local_system() = 0;
52  virtual void fix_velocity(LocalElementAccessorBase<3> ele_ac) = 0;
53  virtual void assemble(LocalElementAccessorBase<3> ele_ac) = 0;
54 
55  // assembly compatible neighbourings
56  virtual void assembly_local_vb(ElementFullIter ele, Neighbour *ngh) = 0;
57 
58  // compute velocity value in the barycenter
59  // TODO: implement and use general interpolations between discrete spaces
61 
63  {}
64 
65 protected:
66 
67  virtual void assemble_sides(LocalElementAccessorBase<3> ele) =0;
68 
70  {}
71 };
72 
73 
74 
75 template <int dim>
77 private:
78  // assembly face integrals (BC)
82 public:
84  : side_quad_(1),
85  fe_p_disc_(0),
86  fe_side_values_(side_map_, side_quad_, fe_p_disc_, update_normal_vectors)
87  {}
89 
90 };
91 
92 
93 
94 template<int dim>
95 class AssemblyMH : public AssemblyBase
96 {
97 public:
99  : quad_(3),
100  fe_values_(map_, quad_, fe_rt_,
102 
103  velocity_interpolation_quad_(0), // veloctiy values in barycenter
104  velocity_interpolation_fv_(map_,velocity_interpolation_quad_, fe_rt_, update_values | update_quadrature_points),
105 
106  ad_(data),
107  loc_system_(size(), size()),
108  loc_system_vb_(2,2)
109  {
110  // local numbering of dofs for MH system
111  unsigned int nsides = dim+1;
112  loc_side_dofs.resize(nsides);
113  loc_ele_dof = nsides;
114  loc_edge_dofs.resize(nsides);
115  for(unsigned int i = 0; i < nsides; i++){
116  loc_side_dofs[i] = i;
117  loc_edge_dofs[i] = nsides + i + 1;
118  }
119  //DebugOut() << print_var(this) << print_var(side_quad_.size());
120 
121  // create local sparsity pattern
122  arma::umat sp(size(),size());
123  sp.zeros();
124  sp.submat(0, 0, nsides, nsides).ones();
125  sp.diag().ones();
126  sp.diag(nsides+1).ones();
127  sp.diag(-nsides-1).ones();
128  loc_system_.set_sparsity(sp);
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);
133  loc_system_vb_.set_sparsity(sp);
134 
135  if (ad_->mortar_method_ == DarcyMH::MortarP0) {
136  mortar_assembly = std::make_shared<P0_CouplingAssembler>(ad_);
137  } else if (ad_->mortar_method_ == DarcyMH::MortarP1) {
138  mortar_assembly = std::make_shared<P1_CouplingAssembler>(ad_);
139  }
140 
141  }
142 
143 
144  ~AssemblyMH<dim>() override
145  {}
146 
147 // LocalSystem& get_local_system() override
148 // { return loc_system_;}
149 
151  {
152  if (mortar_assembly)
153  mortar_assembly->fix_velocity(ele_ac);
154  }
155 
156  void assemble(LocalElementAccessorBase<3> ele_ac) override
157  {
158  ASSERT_EQ_DBG(ele_ac.dim(), dim);
159  loc_system_.reset();
160 
161  set_dofs_and_bc(ele_ac);
162 
163  assemble_sides(ele_ac);
164  assemble_element(ele_ac);
165  assemble_source_term(ele_ac);
166 
167  ad_->lin_sys->set_local_system(loc_system_);
168 
169  assembly_dim_connections(ele_ac);
170 
171  if (ad_->balance != nullptr)
172  add_fluxes_in_balance_matrix(ele_ac);
173 
174  if (mortar_assembly)
175  mortar_assembly->assembly(ele_ac);
176  }
177 
179  {
180  ASSERT_LT_DBG(ele->dim(), 3);
181  //DebugOut() << "alv " << print_var(this);
182  //START_TIMER("Assembly<dim>::assembly_local_vb");
183  // compute normal vector to side
184  arma::vec3 nv;
185  ElementFullIter ele_higher = ad_->mesh->element.full_iter(ngh->side()->element());
186  ngh_values_.fe_side_values_.reinit(ele_higher, ngh->side()->el_idx());
187  nv = ngh_values_.fe_side_values_.normal_vector(0);
188 
189  double value = ad_->sigma.value( ele->centre(), ele->element_accessor()) *
190  2*ad_->conductivity.value( ele->centre(), ele->element_accessor()) *
191  arma::dot(ad_->anisotropy.value( ele->centre(), ele->element_accessor())*nv, nv) *
192  ad_->cross_section.value( ngh->side()->centre(), ele_higher->element_accessor() ) * // cross-section of higher dim. (2d)
193  ad_->cross_section.value( ngh->side()->centre(), ele_higher->element_accessor() ) /
194  ad_->cross_section.value( ele->centre(), ele->element_accessor() ) * // crossection of lower dim.
195  ngh->side()->measure();
196 
197  loc_system_vb_.add_value(0,0, -value);
198  loc_system_vb_.add_value(0,1, value);
199  loc_system_vb_.add_value(1,0, value);
200  loc_system_vb_.add_value(1,1, -value);
201  }
202 
203 
205  {
206  //START_TIMER("Assembly<dim>::make_element_vector");
207  arma::vec3 flux_in_center;
208  flux_in_center.zeros();
209 
210  velocity_interpolation_fv_.reinit(ele);
211  for (unsigned int li = 0; li < ele->n_sides(); li++) {
212  flux_in_center += ad_->mh_dh->side_flux( *(ele->side( li ) ) )
213  * velocity_interpolation_fv_.vector_view(0).value(li,0);
214  }
215 
216  flux_in_center /= ad_->cross_section.value(ele->centre(), ele->element_accessor() );
217  return flux_in_center;
218  }
219 
220 protected:
221  static const unsigned int size()
222  {
223  // dofs: velocity, pressure, edge pressure
225  }
226 
228 
229  ASSERT_DBG(ele_ac.dim() == dim);
230 
231  //set global dof for element (pressure)
232  loc_system_.row_dofs[loc_ele_dof] = loc_system_.col_dofs[loc_ele_dof] = ele_ac.ele_row();
233 
234  //shortcuts
235  const unsigned int nsides = ele_ac.n_sides();
236  LinSys *ls = ad_->lin_sys;
237 
238  Boundary *bcd;
239  unsigned int side_row, edge_row;
240 
241  dirichlet_edge.resize(nsides);
242  for (unsigned int i = 0; i < nsides; i++) {
243 
244  side_row = loc_side_dofs[i]; //local
245  edge_row = loc_edge_dofs[i]; //local
246  loc_system_.row_dofs[side_row] = loc_system_.col_dofs[side_row] = ele_ac.side_row(i); //global
247  loc_system_.row_dofs[edge_row] = loc_system_.col_dofs[edge_row] = ele_ac.edge_row(i); //global
248 
249  bcd = ele_ac.side(i)->cond();
250  dirichlet_edge[i] = 0;
251  if (bcd) {
252  ElementAccessor<3> b_ele = bcd->element_accessor();
253  DarcyMH::EqData::BC_Type type = (DarcyMH::EqData::BC_Type)ad_->bc_type.value(b_ele.centre(), b_ele);
254 
255  double cross_section = ad_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor());
256 
257  if ( type == DarcyMH::EqData::none) {
258  // homogeneous neumann
259  } else if ( type == DarcyMH::EqData::dirichlet ) {
260  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
261  loc_system_.set_solution(loc_edge_dofs[i],bc_pressure,-1);
262  dirichlet_edge[i] = 1;
263 
264  } else if ( type == DarcyMH::EqData::total_flux) {
265  // internally we work with outward flux
266  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
267  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
268  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
269 
270  dirichlet_edge[i] = 2; // to be skipped in LMH source assembly
271  loc_system_.add_value(edge_row, edge_row,
272  -bcd->element()->measure() * bc_sigma * cross_section,
273  (bc_flux - bc_sigma * bc_pressure) * bcd->element()->measure() * cross_section);
274  }
275  else if (type==DarcyMH::EqData::seepage) {
276  ad_->is_linear=false;
277 
278  unsigned int loc_edge_idx = bcd->bc_ele_idx_;
279  char & switch_dirichlet = ad_->bc_switch_dirichlet[loc_edge_idx];
280  double bc_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
281  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
282  double side_flux = bc_flux * bcd->element()->measure() * cross_section;
283 
284  // ** Update BC type. **
285  if (switch_dirichlet) {
286  // check and possibly switch to flux BC
287  // The switch raise error on the corresponding edge row.
288  // Magnitude of the error is abs(solution_flux - side_flux).
289  ASSERT_DBG(ad_->mh_dh->rows_ds->is_local(ele_ac.side_row(i)))(ele_ac.side_row(i));
290  unsigned int loc_side_row = ele_ac.side_local_row(i);
291  double & solution_flux = ls->get_solution_array()[loc_side_row];
292 
293  if ( solution_flux < side_flux) {
294  //DebugOut().fmt("x: {}, to neum, p: {} f: {} -> f: {}\n", b_ele.centre()[0], bc_pressure, solution_flux, side_flux);
295  solution_flux = side_flux;
296  switch_dirichlet=0;
297  }
298  } else {
299  // check and possibly switch to pressure BC
300  // TODO: What is the appropriate DOF in not local?
301  // The switch raise error on the corresponding side row.
302  // Magnitude of the error is abs(solution_head - bc_pressure)
303  // Since usually K is very large, this error would be much
304  // higher then error caused by the inverse switch, this
305  // cause that a solution with the flux violating the
306  // flux inequality leading may be accepted, while the error
307  // in pressure inequality is always satisfied.
308  ASSERT_DBG(ad_->mh_dh->rows_ds->is_local(ele_ac.edge_row(i)))(ele_ac.edge_row(i));
309  unsigned int loc_edge_row = ele_ac.edge_local_row(i);
310  double & solution_head = ls->get_solution_array()[loc_edge_row];
311 
312  if ( solution_head > bc_pressure) {
313  //DebugOut().fmt("x: {}, to dirich, p: {} -> p: {} f: {}\n",b_ele.centre()[0], solution_head, bc_pressure, bc_flux);
314  solution_head = bc_pressure;
315  switch_dirichlet=1;
316  }
317  }
318 
319  // ** Apply BCUpdate BC type. **
320  // Force Dirichlet type during the first iteration of the unsteady case.
321  if (switch_dirichlet || ad_->force_bc_switch ) {
322  //DebugOut().fmt("x: {}, dirich, bcp: {}\n", b_ele.centre()[0], bc_pressure);
323  loc_system_.set_solution(loc_edge_dofs[i],bc_pressure, -1);
324  dirichlet_edge[i] = 1;
325  } else {
326  //DebugOut()("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], side_flux, bc_flux);
327  loc_system_.add_value(edge_row, side_flux);
328  }
329 
330  } else if (type==DarcyMH::EqData::river) {
331  ad_->is_linear=false;
332 
333  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
334  double bc_switch_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
335  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
336  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
337  ASSERT_DBG(ad_->mh_dh->rows_ds->is_local(ele_ac.edge_row(i)))(ele_ac.edge_row(i));
338  unsigned int loc_edge_row = ele_ac.edge_local_row(i);
339  double & solution_head = ls->get_solution_array()[loc_edge_row];
340 
341  // Force Robin type during the first iteration of the unsteady case.
342  if (solution_head > bc_switch_pressure || ad_->force_bc_switch) {
343  // Robin BC
344  //DebugOut().fmt("x: {}, robin, bcp: {}\n", b_ele.centre()[0], bc_pressure);
345  loc_system_.add_value(edge_row, edge_row,
346  -bcd->element()->measure() * bc_sigma * cross_section,
347  bcd->element()->measure() * cross_section * (bc_flux - bc_sigma * bc_pressure) );
348  } else {
349  // Neumann BC
350  //DebugOut().fmt("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], bc_switch_pressure, bc_pressure);
351  double bc_total_flux = bc_flux + bc_sigma*(bc_switch_pressure - bc_pressure);
352 
353  loc_system_.add_value(edge_row, bc_total_flux * bcd->element()->measure() * cross_section);
354  }
355  }
356  else {
357  xprintf(UsrErr, "BC type not supported.\n");
358  }
359  }
360  loc_system_.add_value(side_row, edge_row, 1.0);
361  loc_system_.add_value(edge_row, side_row, 1.0);
362  }
363 
364 // DBGCOUT(<< "ele " << ele_ac.ele_global_idx() << ": ");
365 // for (unsigned int i = 0; i < nsides; i++) cout << loc_system_.row_dofs[loc_side_dofs[i]] << " ";
366 // cout << loc_system_.row_dofs[loc_ele_dof] << " ";
367 // for (unsigned int i = 0; i < nsides; i++) cout << loc_system_.row_dofs[loc_edge_dofs[i]] << " ";
368 // cout << "\n";
369  }
370 
372  {
373  double cs = ad_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor());
374  double conduct = ad_->conductivity.value(ele_ac.centre(), ele_ac.element_accessor());
375  double scale = 1 / cs /conduct;
376 
377  assemble_sides_scale(ele_ac, scale);
378  }
379 
381  {
382  arma::vec3 &gravity_vec = ad_->gravity_vec_;
383 
384  ElementFullIter ele =ele_ac.full_iter();
385  fe_values_.reinit(ele);
386  unsigned int ndofs = fe_values_.get_fe()->n_dofs();
387  unsigned int qsize = fe_values_.get_quadrature()->size();
388  auto velocity = fe_values_.vector_view(0);
389 
390  for (unsigned int k=0; k<qsize; k++)
391  for (unsigned int i=0; i<ndofs; i++){
392  double rhs_val =
393  arma::dot(gravity_vec,velocity.value(i,k))
394  * fe_values_.JxW(k);
395  loc_system_.add_value(i, rhs_val);
396 
397  for (unsigned int j=0; j<ndofs; j++){
398  double mat_val =
399  arma::dot(velocity.value(i,k), //TODO: compute anisotropy before
400  (ad_->anisotropy.value(ele_ac.centre(), ele_ac.element_accessor() )).i()
401  * velocity.value(j,k))
402  * scale * fe_values_.JxW(k);
403 
404  loc_system_.add_value(i, j, mat_val);
405  }
406  }
407 
408  // assemble matrix for weights in BDDCML
409  // approximation to diagonal of
410  // S = -C - B*inv(A)*B'
411  // as
412  // diag(S) ~ - diag(C) - 1./diag(A)
413  // the weights form a partition of unity to average a discontinuous solution from neighbouring subdomains
414  // to a continuous one
415  // it is important to scale the effect - if conductivity is low for one subdomain and high for the other,
416  // trust more the one with low conductivity - it will be closer to the truth than an arithmetic average
417  if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
418  const arma::mat& local_matrix = loc_system_.get_matrix();
419  for(unsigned int i=0; i < ndofs; i++) {
420  double val_side = local_matrix(i,i);
421  double val_edge = -1./local_matrix(i,i);
422 
423  unsigned int side_row = loc_system_.row_dofs[loc_side_dofs[i]];
424  unsigned int edge_row = loc_system_.row_dofs[loc_edge_dofs[i]];
425  static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( side_row, val_side );
426  static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( edge_row, val_edge );
427  }
428  }
429  }
430 
431 
433  // set block B, B': element-side, side-element
434 
435  for(unsigned int side = 0; side < loc_side_dofs.size(); side++){
436  loc_system_.add_value(loc_ele_dof, loc_side_dofs[side], -1.0);
437  loc_system_.add_value(loc_side_dofs[side], loc_ele_dof, -1.0);
438  }
439 
440  if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
441  double val_ele = 1.;
442  static_cast<LinSys_BDDC*>(ad_->lin_sys)->
443  diagonal_weights_set_value( loc_system_.row_dofs[loc_ele_dof], val_ele );
444  }
445  }
446 
447 
449  //D, E',E block: compatible connections: element-edge
450  ElementFullIter ele = ele_ac.full_iter();
451 
452  // no Neighbours => nothing to asssemble here
453  if(ele->n_neighs_vb == 0) return;
454 
455  int ele_row = ele_ac.ele_row();
456  Neighbour *ngh;
457 
458  //DebugOut() << "adc " << print_var(this) << print_var(side_quad_.size());
459  for (unsigned int i = 0; i < ele->n_neighs_vb; i++) {
460  // every compatible connection adds a 2x2 matrix involving
461  // current element pressure and a connected edge pressure
462  ngh = ele_ac.full_iter()->neigh_vb[i];
463  loc_system_vb_.reset();
464  loc_system_vb_.row_dofs[0] = loc_system_vb_.col_dofs[0] = ele_row;
465  loc_system_vb_.row_dofs[1] = loc_system_vb_.col_dofs[1] = ad_->mh_dh->row_4_edge[ ngh->edge_idx() ];
466 
467  assembly_local_vb(ele, ngh);
468 
469  ad_->lin_sys->set_local_system(loc_system_vb_);
470 
471  // update matrix for weights in BDDCML
472  if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
473  int ind = loc_system_vb_.row_dofs[1];
474  // there is -value on diagonal in block C!
475  double new_val = loc_system_vb_.get_matrix()(0,0);
476  static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( ind, new_val );
477  }
478  }
479  }
480 
482 
483  for (unsigned int i = 0; i < ele_ac.n_sides(); i++) {
484  Boundary* bcd = ele_ac.side(i)->cond();
485 
486  if (bcd) {
487  /*
488  DebugOut().fmt("add_flux: {} {} {} {}\n",
489  ad_->mh_dh->el_ds->myp(),
490  ele_ac.ele_global_idx(),
491  ad_->local_boundary_index,
492  ele_ac.side_row(i));
493  */
494  ad_->balance->add_flux_matrix_values(ad_->water_balance_idx, ad_->local_boundary_index,
495  {(IdxInt)(ele_ac.side_row(i))}, {1});
496  ++(ad_->local_boundary_index);
497  }
498  }
499  }
500 
501 
502 
503  // assembly volume integrals
508 
509  NeighSideValues<dim<3?dim:2> ngh_values_;
510 
511  // Interpolation of velocity into barycenters
512  QGauss<dim> velocity_interpolation_quad_;
513  FEValues<dim,3> velocity_interpolation_fv_;
514 
515  // data shared by assemblers of different dimension
516  AssemblyDataPtr ad_;
517  std::vector<unsigned int> dirichlet_edge;
518 
519  LocalSystem loc_system_;
520  LocalSystem loc_system_vb_;
521  std::vector<unsigned int> loc_side_dofs;
522  std::vector<unsigned int> loc_edge_dofs;
523  unsigned int loc_ele_dof;
524 
525  std::shared_ptr<MortarAssemblyBase> mortar_assembly;
526 };
527 
528 
529 #endif /* SRC_FLOW_DARCY_FLOW_ASSEMBLY_HH_ */
const arma::vec3 centre() const
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
double measure() const
Computes the measure of the element.
Definition: elements.cc:90
MappingP1< dim, 3 > map_
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:331
FESideValues< dim+1, 3 > fe_side_values_
std::vector< std::shared_ptr< AssemblyBase > > MultidimAssembly
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void assemble_sides(LocalElementAccessorBase< 3 > ele_ac) override
ElementAccessor< 3 > element_accessor()
virtual arma::vec3 make_element_vector(ElementFullIter ele)=0
void add_fluxes_in_balance_matrix(LocalElementAccessorBase< 3 > ele_ac)
static const unsigned int size()
void assemble(LocalElementAccessorBase< 3 > ele_ac) override
Boundary * cond() const
Definition: side_impl.hh:70
Definition: abscissa.h:25
arma::vec3 make_element_vector(ElementFullIter ele) override
Wrappers for linear systems based on MPIAIJ and MATIS format.
arma::vec3 centre() const
Definition: sides.cc:122
std::shared_ptr< DarcyMH::EqData > AssemblyDataPtr
virtual void assemble_source_term(LocalElementAccessorBase< 3 > ele)
virtual void assemble(LocalElementAccessorBase< 3 > ele_ac)=0
Class FEValues calculates finite element data on the actual cells such as shape function values...
virtual void assembly_local_vb(ElementFullIter ele, Neighbour *ngh)=0
static constexpr bool value
Definition: json.hpp:87
Symmetric Gauss-Legendre quadrature formulae on simplices.
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
std::shared_ptr< DarcyMH::EqData > AssemblyDataPtr
virtual void fix_velocity(LocalElementAccessorBase< 3 > ele_ac)=0
Transformed quadrature points.
FE_RT0< dim, 3 > fe_rt_
double * get_solution_array()
Definition: linsys.hh:305
void set_dofs_and_bc(LocalElementAccessorBase< 3 > ele_ac)
unsigned int edge_idx()
Definitions of basic Lagrangean finite elements with polynomial shape functions.
SideIter side()
void assembly_local_vb(ElementFullIter ele, Neighbour *ngh) override
#define xprintf(...)
Definition: system.hh:92
void assemble_sides_scale(LocalElementAccessorBase< 3 > ele_ac, double scale)
Element * element()
Definition: boundaries.cc:39
Shape function gradients.
Definition: update_flags.hh:95
virtual void update_water_content(LocalElementAccessorBase< 3 > ele)
double measure() const
Definition: sides.cc:29
QGauss< dim > quad_
Normal vectors.
unsigned int bc_ele_idx_
Definition: boundaries.h:76
Definition: system.hh:64
virtual ~AssemblyBase()
QGauss< dim > side_quad_
ElementFullIter element() const
Definition: side_impl.hh:52
FEValues< dim, 3 > fe_values_
void assembly_dim_connections(LocalElementAccessorBase< 3 > ele_ac)
#define ASSERT_DBG(expr)
Definition: asserts.hh:349
Definitions of particular quadrature rules on simplices.
virtual void assemble_sides(LocalElementAccessorBase< 3 > ele)=0
arma::vec::fixed< spacedim > centre() const
Definition: accessors.hh:91
void fix_velocity(LocalElementAccessorBase< 3 > ele_ac) override
unsigned int el_idx() const
Definition: side_impl.hh:81
static MultidimAssembly create(typename Impl< 1 >::AssemblyDataPtr data)
FE_P_disc< dim+1, 3 > fe_p_disc_
void assemble_element(LocalElementAccessorBase< 3 > ele_ac)
MappingP1< dim+1, 3 > side_map_
ElementFullIter full_iter()
Definitions of Raviart-Thomas finite elements.
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:47
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:299