Flow123d  JS_before_hm-915-gc632be9
assembly_lmh.hh
Go to the documentation of this file.
1 /*
2  * assembly_lmh.hh
3  *
4  * Created on: Apr 21, 2016
5  * Author: jb
6  */
7 
8 #ifndef SRC_ASSEMBLY_LMH_HH_
9 #define SRC_ASSEMBLY_LMH_HH_
10 
11 #include "system/index_types.hh"
12 #include "system/fmt/posix.h"
13 #include "mesh/mesh.h"
14 #include "mesh/accessors.hh"
15 #include "mesh/neighbours.h"
16 #include "fem/mapping_p1.hh"
17 #include "fem/fe_p.hh"
18 #include "fem/fe_values.hh"
19 #include "fem/fe_rt.hh"
20 #include "fem/fe_values_views.hh"
21 #include "fem/fe_system.hh"
23 
24 #include "la/linsys_PETSC.hh"
25 // #include "la/linsys_BDDC.hh"
26 #include "la/schur.hh"
27 #include "la/local_system.hh"
28 
29 #include "coupling/balance.hh"
30 #include "flow/assembly_mh.hh"
31 #include "flow/darcy_flow_lmh.hh"
32 #include "flow/mortar_assembly.hh"
33 
34 /** Copy of the assembly class for MH implementation,
35  * with Lumping and further improvements.
36  * Used also for Richards.
37  */
38 template <int dim>
39 class AssemblyLMH : public AssemblyBase
40 {
41 public:
42  typedef std::shared_ptr<DarcyLMH::EqData> AssemblyDataPtrLMH;
43 
44  AssemblyLMH<dim>(AssemblyDataPtrLMH data)
45  : quad_(dim, 3),
46  velocity_interpolation_quad_(dim, 0), // veloctiy values in barycenter
47  ad_(data)
48  {
50  velocity_interpolation_fv_.initialize(velocity_interpolation_quad_, fe_rt_, update_values | update_quadrature_points);
51  // local numbering of dofs for MH system
52  // note: this shortcut supposes that the fe_system is the same on all elements
53  // the function DiscreteSpace.fe(ElementAccessor) does not in fact depend on the element accessor
54  auto fe = ad_->dh_->ds()->fe(ad_->dh_->own_range().begin()->elm())[Dim<dim>{}];
55  FESystem<dim>* fe_system = dynamic_cast<FESystem<dim>*>(fe.get());
56  loc_side_dofs = fe_system->fe_dofs(0);
57  loc_ele_dof = fe_system->fe_dofs(1)[0];
58  loc_edge_dofs = fe_system->fe_dofs(2);
59 
60  FEAL_ASSERT(ad_->mortar_method_ == DarcyFlowInterface::NoMortar)
61  .error("Mortar methods are not supported in Lumped Mixed-Hybrid Method.");
62  // if (ad_->mortar_method_ == DarcyFlowInterface::MortarP0) {
63  // mortar_assembly = std::make_shared<P0_CouplingAssembler>(ad_);
64  // } else if (ad_->mortar_method_ == DarcyFlowInterface::MortarP1) {
65  // mortar_assembly = std::make_shared<P1_CouplingAssembler>(ad_);
66  // }
67 
68  // reconstructed vector for the velocity and pressure
69  // length = local Schur complement offset in LocalSystem
70  schur_offset_ = loc_edge_dofs[0];
71  reconstructed_solution_.zeros(schur_offset_);
72  }
73 
74 
75  ~AssemblyLMH<dim>() override
76  {}
77 
78  void fix_velocity(const DHCellAccessor&) override
79  {
80  // if (mortar_assembly)
81  // mortar_assembly->fix_velocity(ele_ac);
82  }
83 
84  void assemble_reconstruct(const DHCellAccessor& dh_cell) override
85  {
86  ASSERT_EQ_DBG(dh_cell.dim(), dim);
87 
88  assemble_local_system(dh_cell, false); //do not switch dirichlet in seepage when reconstructing
89 
90  // TODO:
91  // if (mortar_assembly)
92  // mortar_assembly->assembly(ele_ac);
93  // if (mortar_assembly)
94  // mortar_assembly->fix_velocity(ele_ac);
95 
96  arma::vec schur_solution = ad_->p_edge_solution.get_subvec(loc_schur_.row_dofs);
97  // reconstruct the velocity and pressure
98  loc_system_.reconstruct_solution_schur(schur_offset_, schur_solution, reconstructed_solution_);
99 
100  // postprocess the velocity
101  postprocess_velocity(dh_cell, reconstructed_solution_);
102 
103  ad_->full_solution.set_subvec(loc_system_.row_dofs.head(schur_offset_), reconstructed_solution_);
104  ad_->full_solution.set_subvec(loc_system_.row_dofs.tail(loc_schur_.row_dofs.n_elem), schur_solution);
105  }
106 
107  void assemble(const DHCellAccessor& dh_cell) override
108  {
109  ASSERT_EQ_DBG(dh_cell.dim(), dim);
110  save_local_system_ = false;
111  bc_fluxes_reconstruted = false;
112 
113  assemble_local_system(dh_cell, true); //do use_dirichlet_switch
114 
115  loc_system_.compute_schur_complement(schur_offset_, loc_schur_, true);
116 
117  save_local_system(dh_cell);
118 
119  loc_schur_.eliminate_solution();
120  ad_->lin_sys_schur->set_local_system(loc_schur_, ad_->dh_cr_->get_local_to_global_map());
121 
122  // TODO:
123  // if (mortar_assembly)
124  // mortar_assembly->assembly(dh_cell);
125  }
126 
127  void update_water_content(const DHCellAccessor&) override
128  {};
129 
130 protected:
131  static unsigned int size()
132  {
133  // dofs: velocity, pressure, edge pressure
135  }
136 
137  void assemble_local_system(const DHCellAccessor& dh_cell, bool use_dirichlet_switch)
138  {
139  set_dofs(dh_cell);
140 
141  assemble_bc(dh_cell, use_dirichlet_switch);
142  assemble_sides(dh_cell);
143  assemble_element(dh_cell);
144  assemble_source_term(dh_cell);
145  assembly_dim_connections(dh_cell);
146  }
147 
148  /** Loads the local system from a map: element index -> LocalSystem,
149  * if it exits, or if the full solution is not yet reconstructed,
150  * and reconstructs the full solution on the element.
151  * Currently used only for seepage BC.
152  */
153  void load_local_system(const DHCellAccessor& dh_cell)
154  {
155  // do this only once per element
156  if(bc_fluxes_reconstruted)
157  return;
158 
159  // possibly find the corresponding local system
160  auto ls = ad_->seepage_bc_systems.find(dh_cell.elm_idx());
161  if (ls != ad_->seepage_bc_systems.end())
162  {
163  arma::vec schur_solution = ad_->p_edge_solution.get_subvec(loc_schur_.row_dofs);
164  // reconstruct the velocity and pressure
165  ls->second.reconstruct_solution_schur(schur_offset_, schur_solution, reconstructed_solution_);
166 
167  postprocess_velocity(dh_cell, reconstructed_solution_);
168 
169  bc_fluxes_reconstruted = true;
170  }
171  }
172 
173 
174  /// Saves the local system to a map: element index -> LocalSystem.
175  /// Currently used only for seepage BC.
176  void save_local_system(const DHCellAccessor& dh_cell)
177  {
178  // for seepage BC, save local system
179  if(save_local_system_)
180  ad_->seepage_bc_systems[dh_cell.elm_idx()] = loc_system_;
181  }
182 
183 
184  void set_dofs(const DHCellAccessor& dh_cell){
185  const ElementAccessor<3> ele = dh_cell.elm();
186  const DHCellAccessor dh_cr_cell = dh_cell.cell_with_other_dh(ad_->dh_cr_.get());
187 
188  unsigned int loc_size = size() + ele->n_neighs_vb();
189  unsigned int loc_size_schur = ele->n_sides() + ele->n_neighs_vb();
190  LocDofVec dofs(loc_size);
191  LocDofVec dofs_schur(loc_size_schur);
192 
193  // add full vec indices
194  dofs.head(dh_cell.n_dofs()) = dh_cell.get_loc_dof_indices();
195  // add schur vec indices
196  dofs_schur.head(dh_cr_cell.n_dofs()) = dh_cr_cell.get_loc_dof_indices();
197 
198  if(ele->n_neighs_vb() != 0)
199  {
200  //D, E',E block: compatible connections: element-edge
201  unsigned int i = 0;
202 
203  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
204 
205  // read neighbor dofs (dh dofhandler)
206  // neighbor cell owning neighb_side
207  DHCellAccessor dh_neighb_cell = neighb_side.cell();
208 
209  // read neighbor dofs (dh_cr dofhandler)
210  // neighbor cell owning neighb_side
211  DHCellAccessor dh_neighb_cr_cell = dh_neighb_cell.cell_with_other_dh(ad_->dh_cr_.get());
212 
213  // local index of pedge dof in local system
214  const unsigned int p = size()+i;
215  // local index of pedge dof on neighboring cell
216  const unsigned int t = dh_neighb_cell.n_dofs() - dh_neighb_cr_cell.n_dofs() + neighb_side.side().side_idx();
217  dofs[p] = dh_neighb_cell.get_loc_dof_indices()[t];
218 
219  // local index of pedge dof in local schur system
220  const unsigned int tt = dh_cr_cell.n_dofs()+i;
221  dofs_schur[tt] = dh_neighb_cr_cell.get_loc_dof_indices()
222  [neighb_side.side().side_idx()];
223  i++;
224  }
225  }
226  loc_system_.reset(dofs, dofs);
227  loc_schur_.reset(dofs_schur, dofs_schur);
228  }
229 
230 
231  void assemble_bc(const DHCellAccessor& dh_cell, bool use_dirichlet_switch){
232  const ElementAccessor<3> ele = dh_cell.elm();
233 
234  dirichlet_edge.resize(ele->n_sides());
235  for(DHCellSide dh_side : dh_cell.side_range()){
236  unsigned int sidx = dh_side.side_idx();
237  dirichlet_edge[sidx] = 0;
238 
239  // assemble BC
240  if (dh_side.side().is_boundary()) {
241  double cross_section = ad_->cross_section.value(ele.centre(), ele);
242  assemble_side_bc(dh_side, cross_section, use_dirichlet_switch);
243 
244  ad_->balance->add_flux_values(ad_->water_balance_idx, dh_side,
245  {loc_system_.row_dofs[loc_side_dofs[sidx]]},
246  {1}, 0);
247  }
248 
249  // side-edge (flux-lambda) terms
250  loc_system_.add_value(loc_side_dofs[sidx], loc_edge_dofs[sidx], 1.0);
251  loc_system_.add_value(loc_edge_dofs[sidx], loc_side_dofs[sidx], 1.0);
252  }
253  }
254 
255 
256  void assemble_side_bc(const DHCellSide& side, double cross_section, bool use_dirichlet_switch)
257  {
258  const unsigned int sidx = side.side_idx();
259  const unsigned int side_row = loc_side_dofs[sidx]; //local
260  const unsigned int edge_row = loc_edge_dofs[sidx]; //local
261 
262  ElementAccessor<3> b_ele = side.cond().element_accessor();
263  DarcyMH::EqData::BC_Type type = (DarcyMH::EqData::BC_Type)ad_->bc_type.value(b_ele.centre(), b_ele);
264 
265  if ( type == DarcyMH::EqData::none) {
266  // homogeneous neumann
267  } else if ( type == DarcyMH::EqData::dirichlet ) {
268  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
269  loc_schur_.set_solution(sidx, bc_pressure);
270  dirichlet_edge[sidx] = 1;
271 
272  } else if ( type == DarcyMH::EqData::total_flux) {
273  // internally we work with outward flux
274  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
275  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
276  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
277 
278  dirichlet_edge[sidx] = 2; // to be skipped in LMH source assembly
279  loc_system_.add_value(edge_row, edge_row,
280  -b_ele.measure() * bc_sigma * cross_section,
281  (bc_flux - bc_sigma * bc_pressure) * b_ele.measure() * cross_section);
282  }
283  else if (type==DarcyMH::EqData::seepage) {
284  ad_->is_linear=false;
285 
286  char & switch_dirichlet = ad_->bc_switch_dirichlet[b_ele.idx()];
287  double bc_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
288  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
289  double side_flux = bc_flux * b_ele.measure() * cross_section;
290 
291  // ** Update BC type. **
292  if(use_dirichlet_switch){ // skip BC change if only reconstructing the solution
293  if (switch_dirichlet) {
294  // check and possibly switch to flux BC
295  // The switch raise error on the corresponding edge row.
296  // Magnitude of the error is abs(solution_flux - side_flux).
297 
298  // try reconstructing local system for seepage BC
299  load_local_system(side.cell());
300  double solution_flux = reconstructed_solution_[side_row];
301 
302  if ( solution_flux < side_flux) {
303  //DebugOut().fmt("x: {}, to neum, p: {} f: {} -> f: {}\n", b_ele.centre()[0], bc_pressure, solution_flux, side_flux);
304  switch_dirichlet=0;
305  }
306  } else {
307  // check and possibly switch to pressure BC
308  // TODO: What is the appropriate DOF in not local?
309  // The switch raise error on the corresponding side row.
310  // Magnitude of the error is abs(solution_head - bc_pressure)
311  // Since usually K is very large, this error would be much
312  // higher then error caused by the inverse switch, this
313  // cause that a solution with the flux violating the
314  // flux inequality leading may be accepted, while the error
315  // in pressure inequality is always satisfied.
316 
317  double solution_head = ad_->p_edge_solution[loc_schur_.row_dofs[sidx]];
318 
319  if ( solution_head > bc_pressure) {
320  //DebugOut().fmt("x: {}, to dirich, p: {} -> p: {} f: {}\n",b_ele.centre()[0], solution_head, bc_pressure, bc_flux);
321  switch_dirichlet=1;
322  }
323  }
324  }
325 
326  save_local_system_ = (bool) switch_dirichlet;
327 
328  // ** Apply BCUpdate BC type. **
329  // Force Dirichlet type during the first iteration of the unsteady case.
330  if (switch_dirichlet || ad_->force_no_neumann_bc ) {
331  //DebugOut().fmt("x: {}, dirich, bcp: {}\n", b_ele.centre()[0], bc_pressure);
332  loc_schur_.set_solution(sidx, bc_pressure);
333  dirichlet_edge[sidx] = 1;
334  } else {
335  //DebugOut()("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], side_flux, bc_flux);
336  loc_system_.add_value(edge_row, side_flux);
337  }
338 
339  } else if (type==DarcyMH::EqData::river) {
340  ad_->is_linear=false;
341 
342  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
343  double bc_switch_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
344  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
345  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
346 
347  double solution_head = ad_->p_edge_solution[loc_schur_.row_dofs[sidx]];
348 
349  // Force Robin type during the first iteration of the unsteady case.
350  if (solution_head > bc_switch_pressure || ad_->force_no_neumann_bc) {
351  // Robin BC
352  //DebugOut().fmt("x: {}, robin, bcp: {}\n", b_ele.centre()[0], bc_pressure);
353  loc_system_.add_value(edge_row, edge_row,
354  -b_ele.measure() * bc_sigma * cross_section,
355  b_ele.measure() * cross_section * (bc_flux - bc_sigma * bc_pressure) );
356  } else {
357  // Neumann BC
358  //DebugOut().fmt("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], bc_switch_pressure, bc_pressure);
359  double bc_total_flux = bc_flux + bc_sigma*(bc_switch_pressure - bc_pressure);
360 
361  loc_system_.add_value(edge_row, bc_total_flux * b_ele.measure() * cross_section);
362  }
363  }
364  else {
365  xprintf(UsrErr, "BC type not supported.\n");
366  }
367  }
368 
369 
370  virtual void assemble_sides(const DHCellAccessor& dh_cell)
371  {
372  const ElementAccessor<3> ele = dh_cell.elm();
373  double cs = ad_->cross_section.value(ele.centre(), ele);
374  double conduct = ad_->conductivity.value(ele.centre(), ele);
375  double scale = 1 / cs /conduct;
376 
377  assemble_sides_scale(dh_cell, scale);
378  }
379 
380  void assemble_sides_scale(const DHCellAccessor& dh_cell, double scale)
381  {
382  arma::vec3 &gravity_vec = ad_->gravity_vec_;
383  auto ele = dh_cell.elm();
384 
385  fe_values_.reinit(ele);
386  unsigned int ndofs = fe_values_.n_dofs();
387  unsigned int qsize = fe_values_.n_points();
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.centre(), ele)).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  virtual void assemble_source_term(const DHCellAccessor& dh_cell)
448  {
449  const ElementAccessor<3> ele = dh_cell.elm();
450 
451  // compute lumped source
452  double alpha = 1.0 / ele->n_sides();
453  double cross_section = ad_->cross_section.value(ele.centre(), ele);
454  double coef = alpha * ele.measure() * cross_section;
455 
456  double source = ad_->water_source_density.value(ele.centre(), ele)
457  + ad_->extra_source.value(ele.centre(), ele);
458  double source_term = coef * source;
459 
460  // in unsteady, compute time term
461  double storativity = 0.0;
462  double time_term_diag = 0.0, time_term = 0.0, time_term_rhs = 0.0;
463 
464  if(! ad_->use_steady_assembly_)
465  {
466  storativity = ad_->storativity.value(ele.centre(), ele)
467  + ad_->extra_storativity.value(ele.centre(), ele);
468  time_term = coef * storativity;
469  }
470 
471  for (unsigned int i=0; i<ele->n_sides(); i++)
472  {
473  if(! ad_->use_steady_assembly_)
474  {
475  time_term_diag = time_term / ad_->time_step_;
476  time_term_rhs = time_term_diag * ad_->p_edge_solution_previous_time[loc_schur_.row_dofs[i]];
477 
478  ad_->balance->add_mass_values(ad_->water_balance_idx, dh_cell,
479  {loc_system_.row_dofs[loc_edge_dofs[i]]}, {time_term}, 0);
480  }
481 
482  this->loc_system_.add_value(loc_edge_dofs[i], loc_edge_dofs[i],
483  -time_term_diag,
484  -source_term - time_term_rhs);
485 
486  ad_->balance->add_source_values(ad_->water_balance_idx, ele.region().bulk_idx(),
487  {loc_system_.row_dofs[loc_edge_dofs[i]]}, {0},{source_term});
488  }
489  }
490 
492  //D, E',E block: compatible connections: element-edge
493  const ElementAccessor<3> ele = dh_cell.elm();
494 
495  // no Neighbours => nothing to asssemble here
496  if(ele->n_neighs_vb() == 0) return;
497 
498  ASSERT_LT_DBG(ele->dim(), 3);
499  arma::vec3 nv;
500 
501  unsigned int i = 0;
502  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
503  // every compatible connection adds a 2x2 matrix involving
504  // current element pressure and a connected edge pressure
505  unsigned int p = size()+i; // loc dof of higher ele edge
506 
507  ElementAccessor<3> ele_higher = neighb_side.cell().elm();
508  ngh_values_.fe_side_values_.reinit(neighb_side.side());
509  nv = ngh_values_.fe_side_values_.normal_vector(0);
510 
511  double value = ad_->sigma.value( ele.centre(), ele) *
512  2*ad_->conductivity.value( ele.centre(), ele) *
513  arma::dot(ad_->anisotropy.value( ele.centre(), ele)*nv, nv) *
514  ad_->cross_section.value( neighb_side.centre(), ele_higher ) * // cross-section of higher dim. (2d)
515  ad_->cross_section.value( neighb_side.centre(), ele_higher ) /
516  ad_->cross_section.value( ele.centre(), ele ) * // crossection of lower dim.
517  neighb_side.measure();
518 
519  loc_system_.add_value(loc_ele_dof, loc_ele_dof, -value);
520  loc_system_.add_value(loc_ele_dof, p, value);
521  loc_system_.add_value(p,loc_ele_dof, value);
522  loc_system_.add_value(p,p, -value);
523 
524 // // update matrix for weights in BDDCML
525 // if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
526 // int ind = loc_system_.row_dofs[p];
527 // // there is -value on diagonal in block C!
528 // static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( ind, -value );
529 // }
530  ++i;
531  }
532  }
533 
534 
535  void postprocess_velocity(const DHCellAccessor& dh_cell, arma::vec& solution)
536  {
537  const ElementAccessor<3> ele = dh_cell.elm();
538 
539  double edge_scale = ele.measure()
540  * ad_->cross_section.value(ele.centre(), ele)
541  / ele->n_sides();
542 
543  double edge_source_term = edge_scale *
544  ( ad_->water_source_density.value(ele.centre(), ele)
545  + ad_->extra_source.value(ele.centre(), ele));
546 
547  postprocess_velocity_specific(dh_cell, solution, edge_scale, edge_source_term);
548  }
549 
550  virtual void postprocess_velocity_specific(const DHCellAccessor& dh_cell, arma::vec& solution,
551  double edge_scale, double edge_source_term)// override
552  {
553  const ElementAccessor<3> ele = dh_cell.elm();
554 
555  double storativity = ad_->storativity.value(ele.centre(), ele)
556  + ad_->extra_storativity.value(ele.centre(), ele);
557  double new_pressure, old_pressure, time_term = 0.0;
558 
559  for (unsigned int i=0; i<ele->n_sides(); i++) {
560 
561  if( ! ad_->use_steady_assembly_)
562  {
563  new_pressure = ad_->p_edge_solution[loc_schur_.row_dofs[i]];
564  old_pressure = ad_->p_edge_solution_previous_time[loc_schur_.row_dofs[i]];
565  time_term = edge_scale * storativity / ad_->time_step_ * (new_pressure - old_pressure);
566  }
567  solution[loc_side_dofs[i]] += edge_source_term - time_term;
568  }
569  }
570 
571  // assembly volume integrals
575 
576  NeighSideValues<dim<3?dim:2> ngh_values_;
577 
578  // Interpolation of velocity into barycenters
579  QGauss velocity_interpolation_quad_;
580  FEValues<3> velocity_interpolation_fv_;
581 
582  // data shared by assemblers of different dimension
583  AssemblyDataPtrLMH ad_;
584 
585  /** TODO: Investigate why the hell do we need this flag.
586  * If removed, it does not break any of the integration tests,
587  * however it must influence the Dirichlet rows in matrix.
588  */
589  std::vector<unsigned int> dirichlet_edge;
590 
591  LocalSystem loc_system_;
592  LocalSystem loc_schur_;
593  std::vector<unsigned int> loc_side_dofs;
594  std::vector<unsigned int> loc_edge_dofs;
595  unsigned int loc_ele_dof;
596 
597  // std::shared_ptr<MortarAssemblyBase> mortar_assembly;
598 
599  /// Index offset in the local system for the Schur complement.
600  unsigned int schur_offset_;
601 
602  /// Vector for reconstruted solution (velocity and pressure on element) from Schur complement.
603  arma::vec reconstructed_solution_;
604 
605  /// Flag for saving the local system.
606  /// Currently used only in case of seepage BC.
607  bool save_local_system_;
608 
609  /// Flag indicating whether the fluxes for seepage BC has been reconstructed already.
610  bool bc_fluxes_reconstruted;
611 };
612 
613 
614 #endif /* SRC_ASSEMBLY_LMH_HH_ */
virtual void assemble_sides(const DHCellAccessor &dh_cell)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
void assemble_local_system(const DHCellAccessor &dh_cell, bool use_dirichlet_switch)
Definition: mixed.hh:25
Range< DHCellSide > side_range() const
Returns range of cell sides.
static unsigned int size()
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
Definition: fe_values.hh:254
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
ArmaVec< double, N > vec
Definition: armor.hh:861
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void save_local_system(const DHCellAccessor &dh_cell)
FEValues< 3 > fe_values_
unsigned int dim() const
Return dimension of element appropriate to cell.
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
unsigned int n_dofs() const
Returns the number of shape functions.
Definition: fe_values.hh:279
Cell accessor allow iterate over DOF handler cells.
Class FEValues calculates finite element data on the actual cells such as shape function values...
Class FESystem for compound finite elements.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
unsigned int side_idx() const
unsigned int n_dofs() const
Return number of dofs on given cell.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:553
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...
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
void load_local_system(const DHCellAccessor &dh_cell)
unsigned int dim() const
Definition: elements.h:121
std::shared_ptr< DarcyLMH::EqData > AssemblyDataPtrLMH
Definition: assembly_lmh.hh:42
virtual void assemble_source_term(const DHCellAccessor &dh_cell)
Transformed quadrature points.
void assemble_side_bc(const DHCellSide &side, double cross_section, bool use_dirichlet_switch)
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
#define FMT_UNUSED
Definition: posix.h:75
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
void assemble(const DHCellAccessor &dh_cell) override
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:60
Definitions of basic Lagrangean finite elements with polynomial shape functions.
std::vector< unsigned int > fe_dofs(unsigned int fe_index)
Return dof indices belonging to given sub-FE.
Definition: fe_system.cc:267
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...
Boundary cond() const
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
Definition: fe_values.cc:505
void assembly_dim_connections(const DHCellAccessor &dh_cell)
void assemble_element(FMT_UNUSED const DHCellAccessor &dh_cell)
double measure() const
Computes the measure of the element.
Shape function gradients.
Definition: update_flags.hh:95
Region region() const
Definition: accessors.hh:165
FE_RT0< dim > fe_rt_
Definition: system.hh:65
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
virtual void postprocess_velocity_specific(const DHCellAccessor &dh_cell, arma::vec &solution, double edge_scale, double edge_source_term)
void fix_velocity(const DHCellAccessor &) override
Definition: assembly_lmh.hh:78
void assemble_bc(const DHCellAccessor &dh_cell, bool use_dirichlet_switch)
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
Definitions of particular quadrature rules on simplices.
Calculates finite element data on the actual cell.
Definition: fe_values.hh:59
void set_dofs(const DHCellAccessor &dh_cell)
DHCellAccessor cell_with_other_dh(const DOFHandlerMultiDim *dh) const
Create new accessor with same local idx and given DOF handler. Actual and given DOF handler must be c...
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:200
unsigned int n_neighs_vb() const
Return number of neighbours.
Definition: elements.h:67
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
void postprocess_velocity(const DHCellAccessor &dh_cell, arma::vec &solution)
void assemble_reconstruct(const DHCellAccessor &dh_cell) override
Definition: assembly_lmh.hh:84
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:273
Definitions of Raviart-Thomas finite elements.
Side accessor allows to iterate over sides of DOF handler cell.
void update_water_content(const DHCellAccessor &) override
Updates water content in Richards.
ElementAccessor< 3 > element_accessor()
#define FEAL_ASSERT(expr)
Definition of assert for debug and release mode.
Definition: asserts.hh:280
Transformed quadrature weights.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
void assemble_sides_scale(const DHCellAccessor &dh_cell, double scale)