Flow123d  JS_before_hm-1626-gde32303
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  // TODO the function should be DiscreteSpace.fe(ElementAccessor)
54  // and correct form fe(ad_->dh_->own_range().begin()->elm()) (see documentation of DiscreteSpace::fe)
55  auto fe = ad_->dh_->ds()->fe()[Dim<dim>{}];
56  FESystem<dim>* fe_system = dynamic_cast<FESystem<dim>*>(fe.get());
57  loc_side_dofs = fe_system->fe_dofs(0);
58  loc_ele_dof = fe_system->fe_dofs(1)[0];
59  loc_edge_dofs = fe_system->fe_dofs(2);
60 
61  FEAL_ASSERT(ad_->mortar_method_ == DarcyFlowInterface::NoMortar)
62  .error("Mortar methods are not supported in Lumped Mixed-Hybrid Method.");
63  // if (ad_->mortar_method_ == DarcyFlowInterface::MortarP0) {
64  // mortar_assembly = std::make_shared<P0_CouplingAssembler>(ad_);
65  // } else if (ad_->mortar_method_ == DarcyFlowInterface::MortarP1) {
66  // mortar_assembly = std::make_shared<P1_CouplingAssembler>(ad_);
67  // }
68 
69  // reconstructed vector for the velocity and pressure
70  // length = local Schur complement offset in LocalSystem
71  schur_offset_ = loc_edge_dofs[0];
72  reconstructed_solution_.zeros(schur_offset_);
73  }
74 
75 
76  ~AssemblyLMH<dim>() override
77  {}
78 
79  void fix_velocity(const DHCellAccessor&) override
80  {
81  // if (mortar_assembly)
82  // mortar_assembly->fix_velocity(ele_ac);
83  }
84 
85  void assemble_reconstruct(const DHCellAccessor& dh_cell) override
86  {
87  ASSERT_EQ_DBG(dh_cell.dim(), dim);
88 
89  assemble_local_system(dh_cell, false); //do not switch dirichlet in seepage when reconstructing
90 
91  // TODO:
92  // if (mortar_assembly)
93  // mortar_assembly->assembly(ele_ac);
94  // if (mortar_assembly)
95  // mortar_assembly->fix_velocity(ele_ac);
96 
97  arma::vec schur_solution = ad_->p_edge_solution.get_subvec(loc_schur_.row_dofs);
98  // reconstruct the velocity and pressure
99  loc_system_.reconstruct_solution_schur(schur_offset_, schur_solution, reconstructed_solution_);
100 
101  // postprocess the velocity
102  postprocess_velocity(dh_cell, reconstructed_solution_);
103 
104  ad_->full_solution.set_subvec(loc_system_.row_dofs.head(schur_offset_), reconstructed_solution_);
105  ad_->full_solution.set_subvec(loc_system_.row_dofs.tail(loc_schur_.row_dofs.n_elem), schur_solution);
106  }
107 
108  void assemble(const DHCellAccessor& dh_cell) override
109  {
110  ASSERT_EQ_DBG(dh_cell.dim(), dim);
111  save_local_system_ = false;
112  bc_fluxes_reconstruted = false;
113 
114  assemble_local_system(dh_cell, true); //do use_dirichlet_switch
115 
116  loc_system_.compute_schur_complement(schur_offset_, loc_schur_, true);
117 
118  save_local_system(dh_cell);
119 
120  loc_schur_.eliminate_solution();
121  ad_->lin_sys_schur->set_local_system(loc_schur_, ad_->dh_cr_->get_local_to_global_map());
122 
123  // TODO:
124  // if (mortar_assembly)
125  // mortar_assembly->assembly(dh_cell);
126  }
127 
128  void update_water_content(const DHCellAccessor&) override
129  {};
130 
131 protected:
132  static unsigned int size()
133  {
134  // dofs: velocity, pressure, edge pressure
136  }
137 
138  void assemble_local_system(const DHCellAccessor& dh_cell, bool use_dirichlet_switch)
139  {
140  set_dofs(dh_cell);
141 
142  assemble_bc(dh_cell, use_dirichlet_switch);
143  assemble_sides(dh_cell);
144  assemble_element(dh_cell);
145  assemble_source_term(dh_cell);
146  assembly_dim_connections(dh_cell);
147  }
148 
149  /** Loads the local system from a map: element index -> LocalSystem,
150  * if it exits, or if the full solution is not yet reconstructed,
151  * and reconstructs the full solution on the element.
152  * Currently used only for seepage BC.
153  */
154  void load_local_system(const DHCellAccessor& dh_cell)
155  {
156  // do this only once per element
157  if(bc_fluxes_reconstruted)
158  return;
159 
160  // possibly find the corresponding local system
161  auto ls = ad_->seepage_bc_systems.find(dh_cell.elm_idx());
162  if (ls != ad_->seepage_bc_systems.end())
163  {
164  arma::vec schur_solution = ad_->p_edge_solution.get_subvec(loc_schur_.row_dofs);
165  // reconstruct the velocity and pressure
166  ls->second.reconstruct_solution_schur(schur_offset_, schur_solution, reconstructed_solution_);
167 
168  postprocess_velocity(dh_cell, reconstructed_solution_);
169 
170  bc_fluxes_reconstruted = true;
171  }
172  }
173 
174 
175  /// Saves the local system to a map: element index -> LocalSystem.
176  /// Currently used only for seepage BC.
177  void save_local_system(const DHCellAccessor& dh_cell)
178  {
179  // for seepage BC, save local system
180  if(save_local_system_)
181  ad_->seepage_bc_systems[dh_cell.elm_idx()] = loc_system_;
182  }
183 
184 
185  void set_dofs(const DHCellAccessor& dh_cell){
186  const ElementAccessor<3> ele = dh_cell.elm();
187  const DHCellAccessor dh_cr_cell = dh_cell.cell_with_other_dh(ad_->dh_cr_.get());
188 
189  unsigned int loc_size = size() + ele->n_neighs_vb();
190  unsigned int loc_size_schur = ele->n_sides() + ele->n_neighs_vb();
191  LocDofVec dofs(loc_size);
192  LocDofVec dofs_schur(loc_size_schur);
193 
194  // add full vec indices
195  dofs.head(dh_cell.n_dofs()) = dh_cell.get_loc_dof_indices();
196  // add schur vec indices
197  dofs_schur.head(dh_cr_cell.n_dofs()) = dh_cr_cell.get_loc_dof_indices();
198 
199  if(ele->n_neighs_vb() != 0)
200  {
201  //D, E',E block: compatible connections: element-edge
202  unsigned int i = 0;
203 
204  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
205 
206  // read neighbor dofs (dh dofhandler)
207  // neighbor cell owning neighb_side
208  DHCellAccessor dh_neighb_cell = neighb_side.cell();
209 
210  // read neighbor dofs (dh_cr dofhandler)
211  // neighbor cell owning neighb_side
212  DHCellAccessor dh_neighb_cr_cell = dh_neighb_cell.cell_with_other_dh(ad_->dh_cr_.get());
213 
214  // local index of pedge dof in local system
215  const unsigned int p = size()+i;
216  // local index of pedge dof on neighboring cell
217  const unsigned int t = dh_neighb_cell.n_dofs() - dh_neighb_cr_cell.n_dofs() + neighb_side.side().side_idx();
218  dofs[p] = dh_neighb_cell.get_loc_dof_indices()[t];
219 
220  // local index of pedge dof in local schur system
221  const unsigned int tt = dh_cr_cell.n_dofs()+i;
222  dofs_schur[tt] = dh_neighb_cr_cell.get_loc_dof_indices()
223  [neighb_side.side().side_idx()];
224  i++;
225  }
226  }
227  loc_system_.reset(dofs, dofs);
228  loc_schur_.reset(dofs_schur, dofs_schur);
229  }
230 
231 
232  void assemble_bc(const DHCellAccessor& dh_cell, bool use_dirichlet_switch){
233  const ElementAccessor<3> ele = dh_cell.elm();
234 
235  dirichlet_edge.resize(ele->n_sides());
236  for(DHCellSide dh_side : dh_cell.side_range()){
237  unsigned int sidx = dh_side.side_idx();
238  dirichlet_edge[sidx] = 0;
239 
240  // assemble BC
241  if (dh_side.side().is_boundary()) {
242  double cross_section = ad_->cross_section.value(ele.centre(), ele);
243  assemble_side_bc(dh_side, cross_section, use_dirichlet_switch);
244 
245  ad_->balance->add_flux_values(ad_->water_balance_idx, dh_side,
246  {loc_system_.row_dofs[loc_side_dofs[sidx]]},
247  {1}, 0);
248  }
249 
250  // side-edge (flux-lambda) terms
251  loc_system_.add_value(loc_side_dofs[sidx], loc_edge_dofs[sidx], 1.0);
252  loc_system_.add_value(loc_edge_dofs[sidx], loc_side_dofs[sidx], 1.0);
253  }
254  }
255 
256 
257  void assemble_side_bc(const DHCellSide& side, double cross_section, bool use_dirichlet_switch)
258  {
259  const unsigned int sidx = side.side_idx();
260  const unsigned int side_row = loc_side_dofs[sidx]; //local
261  const unsigned int edge_row = loc_edge_dofs[sidx]; //local
262 
263  ElementAccessor<3> b_ele = side.cond().element_accessor();
264  DarcyMH::EqData::BC_Type type = (DarcyMH::EqData::BC_Type)ad_->bc_type.value(b_ele.centre(), b_ele);
265 
266  if ( type == DarcyMH::EqData::none) {
267  // homogeneous neumann
268  } else if ( type == DarcyMH::EqData::dirichlet ) {
269  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
270  loc_schur_.set_solution(sidx, bc_pressure);
271  dirichlet_edge[sidx] = 1;
272 
273  } else if ( type == DarcyMH::EqData::total_flux) {
274  // internally we work with outward flux
275  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
276  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
277  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
278 
279  dirichlet_edge[sidx] = 2; // to be skipped in LMH source assembly
280  loc_system_.add_value(edge_row, edge_row,
281  -b_ele.measure() * bc_sigma * cross_section,
282  (bc_flux - bc_sigma * bc_pressure) * b_ele.measure() * cross_section);
283  }
284  else if (type==DarcyMH::EqData::seepage) {
285  ad_->is_linear=false;
286 
287  char & switch_dirichlet = ad_->bc_switch_dirichlet[b_ele.idx()];
288  double bc_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
289  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
290  double side_flux = bc_flux * b_ele.measure() * cross_section;
291 
292  // ** Update BC type. **
293  if(use_dirichlet_switch){ // skip BC change if only reconstructing the solution
294  if (switch_dirichlet) {
295  // check and possibly switch to flux BC
296  // The switch raise error on the corresponding edge row.
297  // Magnitude of the error is abs(solution_flux - side_flux).
298 
299  // try reconstructing local system for seepage BC
300  load_local_system(side.cell());
301  double solution_flux = reconstructed_solution_[side_row];
302 
303  if ( solution_flux < side_flux) {
304  //DebugOut().fmt("x: {}, to neum, p: {} f: {} -> f: {}\n", b_ele.centre()[0], bc_pressure, solution_flux, side_flux);
305  switch_dirichlet=0;
306  }
307  } else {
308  // check and possibly switch to pressure BC
309  // TODO: What is the appropriate DOF in not local?
310  // The switch raise error on the corresponding side row.
311  // Magnitude of the error is abs(solution_head - bc_pressure)
312  // Since usually K is very large, this error would be much
313  // higher then error caused by the inverse switch, this
314  // cause that a solution with the flux violating the
315  // flux inequality leading may be accepted, while the error
316  // in pressure inequality is always satisfied.
317 
318  double solution_head = ad_->p_edge_solution.get(loc_schur_.row_dofs[sidx]);
319 
320  if ( solution_head > bc_pressure) {
321  //DebugOut().fmt("x: {}, to dirich, p: {} -> p: {} f: {}\n",b_ele.centre()[0], solution_head, bc_pressure, bc_flux);
322  switch_dirichlet=1;
323  }
324  }
325  }
326 
327  save_local_system_ = (bool) switch_dirichlet;
328 
329  // ** Apply BCUpdate BC type. **
330  // Force Dirichlet type during the first iteration of the unsteady case.
331  if (switch_dirichlet || ad_->force_no_neumann_bc ) {
332  //DebugOut().fmt("x: {}, dirich, bcp: {}\n", b_ele.centre()[0], bc_pressure);
333  loc_schur_.set_solution(sidx, bc_pressure);
334  dirichlet_edge[sidx] = 1;
335  } else {
336  //DebugOut()("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], side_flux, bc_flux);
337  loc_system_.add_value(edge_row, side_flux);
338  }
339 
340  } else if (type==DarcyMH::EqData::river) {
341  ad_->is_linear=false;
342 
343  double bc_pressure = ad_->bc_pressure.value(b_ele.centre(), b_ele);
344  double bc_switch_pressure = ad_->bc_switch_pressure.value(b_ele.centre(), b_ele);
345  double bc_flux = -ad_->bc_flux.value(b_ele.centre(), b_ele);
346  double bc_sigma = ad_->bc_robin_sigma.value(b_ele.centre(), b_ele);
347 
348  double solution_head = ad_->p_edge_solution.get(loc_schur_.row_dofs[sidx]);
349 
350  // Force Robin type during the first iteration of the unsteady case.
351  if (solution_head > bc_switch_pressure || ad_->force_no_neumann_bc) {
352  // Robin BC
353  //DebugOut().fmt("x: {}, robin, bcp: {}\n", b_ele.centre()[0], bc_pressure);
354  loc_system_.add_value(edge_row, edge_row,
355  -b_ele.measure() * bc_sigma * cross_section,
356  b_ele.measure() * cross_section * (bc_flux - bc_sigma * bc_pressure) );
357  } else {
358  // Neumann BC
359  //DebugOut().fmt("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], bc_switch_pressure, bc_pressure);
360  double bc_total_flux = bc_flux + bc_sigma*(bc_switch_pressure - bc_pressure);
361 
362  loc_system_.add_value(edge_row, bc_total_flux * b_ele.measure() * cross_section);
363  }
364  }
365  else {
366  THROW( ExcBCNotSupported() );
367  }
368  }
369 
370 
371  virtual void assemble_sides(const DHCellAccessor& dh_cell)
372  {
373  const ElementAccessor<3> ele = dh_cell.elm();
374  double cs = ad_->cross_section.value(ele.centre(), ele);
375  double conduct = ad_->conductivity.value(ele.centre(), ele);
376  double scale = 1 / cs /conduct;
377 
378  assemble_sides_scale(dh_cell, scale);
379  }
380 
381  void assemble_sides_scale(const DHCellAccessor& dh_cell, double scale)
382  {
383  arma::vec3 &gravity_vec = ad_->gravity_vec_;
384  auto ele = dh_cell.elm();
385 
386  fe_values_.reinit(ele);
387  unsigned int ndofs = fe_values_.n_dofs();
388  unsigned int qsize = fe_values_.n_points();
389  auto velocity = fe_values_.vector_view(0);
390 
391  for (unsigned int k=0; k<qsize; k++)
392  for (unsigned int i=0; i<ndofs; i++){
393  double rhs_val =
394  arma::dot(gravity_vec,velocity.value(i,k))
395  * fe_values_.JxW(k);
396  loc_system_.add_value(i, rhs_val);
397 
398  for (unsigned int j=0; j<ndofs; j++){
399  double mat_val =
400  arma::dot(velocity.value(i,k), //TODO: compute anisotropy before
401  (ad_->anisotropy.value(ele.centre(), ele)).i()
402  * velocity.value(j,k))
403  * scale * fe_values_.JxW(k);
404 
405  loc_system_.add_value(i, j, mat_val);
406  }
407  }
408 
409  // assemble matrix for weights in BDDCML
410  // approximation to diagonal of
411  // S = -C - B*inv(A)*B'
412  // as
413  // diag(S) ~ - diag(C) - 1./diag(A)
414  // the weights form a partition of unity to average a discontinuous solution from neighbouring subdomains
415  // to a continuous one
416  // it is important to scale the effect - if conductivity is low for one subdomain and high for the other,
417  // trust more the one with low conductivity - it will be closer to the truth than an arithmetic average
418 // if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
419 // const arma::mat& local_matrix = loc_system_.get_matrix();
420 // for(unsigned int i=0; i < ndofs; i++) {
421 // double val_side = local_matrix(i,i);
422 // double val_edge = -1./local_matrix(i,i);
423 //
424 // unsigned int side_row = loc_system_.row_dofs[loc_side_dofs[i]];
425 // unsigned int edge_row = loc_system_.row_dofs[loc_edge_dofs[i]];
426 // static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( side_row, val_side );
427 // static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( edge_row, val_edge );
428 // }
429 // }
430  }
431 
432 
434  // set block B, B': element-side, side-element
435 
436  for(unsigned int side = 0; side < loc_side_dofs.size(); side++){
437  loc_system_.add_value(loc_ele_dof, loc_side_dofs[side], -1.0);
438  loc_system_.add_value(loc_side_dofs[side], loc_ele_dof, -1.0);
439  }
440 
441 // if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
442 // double val_ele = 1.;
443 // static_cast<LinSys_BDDC*>(ad_->lin_sys)->
444 // diagonal_weights_set_value( loc_system_.row_dofs[loc_ele_dof], val_ele );
445 // }
446  }
447 
448  virtual void assemble_source_term(const DHCellAccessor& dh_cell)
449  {
450  const ElementAccessor<3> ele = dh_cell.elm();
451 
452  // compute lumped source
453  double alpha = 1.0 / ele->n_sides();
454  double cross_section = ad_->cross_section.value(ele.centre(), ele);
455  double coef = alpha * ele.measure() * cross_section;
456 
457  double source = ad_->water_source_density.value(ele.centre(), ele)
458  + ad_->extra_source.value(ele.centre(), ele);
459  double source_term = coef * source;
460 
461  // in unsteady, compute time term
462  double storativity = 0.0;
463  double time_term_diag = 0.0, time_term = 0.0, time_term_rhs = 0.0;
464 
465  if(! ad_->use_steady_assembly_)
466  {
467  storativity = ad_->storativity.value(ele.centre(), ele)
468  + ad_->extra_storativity.value(ele.centre(), ele);
469  time_term = coef * storativity;
470  }
471 
472  for (unsigned int i=0; i<ele->n_sides(); i++)
473  {
474  if(! ad_->use_steady_assembly_)
475  {
476  time_term_diag = time_term / ad_->time_step_;
477  time_term_rhs = time_term_diag * ad_->p_edge_solution_previous_time.get(loc_schur_.row_dofs[i]);
478 
479  ad_->balance->add_mass_values(ad_->water_balance_idx, dh_cell,
480  {loc_system_.row_dofs[loc_edge_dofs[i]]}, {time_term}, 0);
481  }
482 
483  this->loc_system_.add_value(loc_edge_dofs[i], loc_edge_dofs[i],
484  -time_term_diag,
485  -source_term - time_term_rhs);
486 
487  ad_->balance->add_source_values(ad_->water_balance_idx, ele.region().bulk_idx(),
488  {loc_system_.row_dofs[loc_edge_dofs[i]]}, {0},{source_term});
489  }
490  }
491 
493  //D, E',E block: compatible connections: element-edge
494  const ElementAccessor<3> ele = dh_cell.elm();
495 
496  // no Neighbours => nothing to asssemble here
497  if(ele->n_neighs_vb() == 0) return;
498 
499  ASSERT_LT_DBG(ele->dim(), 3);
500  arma::vec3 nv;
501 
502  unsigned int i = 0;
503  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
504  // every compatible connection adds a 2x2 matrix involving
505  // current element pressure and a connected edge pressure
506  unsigned int p = size()+i; // loc dof of higher ele edge
507 
508  ElementAccessor<3> ele_higher = neighb_side.cell().elm();
509  ngh_values_.fe_side_values_.reinit(neighb_side.side());
510  nv = ngh_values_.fe_side_values_.normal_vector(0);
511 
512  double value = ad_->sigma.value( ele.centre(), ele) *
513  2*ad_->conductivity.value( ele.centre(), ele) *
514  arma::dot(ad_->anisotropy.value( ele.centre(), ele)*nv, nv) *
515  ad_->cross_section.value( neighb_side.centre(), ele_higher ) * // cross-section of higher dim. (2d)
516  ad_->cross_section.value( neighb_side.centre(), ele_higher ) /
517  ad_->cross_section.value( ele.centre(), ele ) * // crossection of lower dim.
518  neighb_side.measure();
519 
520  loc_system_.add_value(loc_ele_dof, loc_ele_dof, -value);
521  loc_system_.add_value(loc_ele_dof, p, value);
522  loc_system_.add_value(p,loc_ele_dof, value);
523  loc_system_.add_value(p,p, -value);
524 
525 // // update matrix for weights in BDDCML
526 // if ( typeid(*ad_->lin_sys) == typeid(LinSys_BDDC) ) {
527 // int ind = loc_system_.row_dofs[p];
528 // // there is -value on diagonal in block C!
529 // static_cast<LinSys_BDDC*>(ad_->lin_sys)->diagonal_weights_set_value( ind, -value );
530 // }
531  ++i;
532  }
533  }
534 
535 
536  void postprocess_velocity(const DHCellAccessor& dh_cell, arma::vec& solution)
537  {
538  const ElementAccessor<3> ele = dh_cell.elm();
539 
540  double edge_scale = ele.measure()
541  * ad_->cross_section.value(ele.centre(), ele)
542  / ele->n_sides();
543 
544  double edge_source_term = edge_scale *
545  ( ad_->water_source_density.value(ele.centre(), ele)
546  + ad_->extra_source.value(ele.centre(), ele));
547 
548  postprocess_velocity_specific(dh_cell, solution, edge_scale, edge_source_term);
549  }
550 
551  virtual void postprocess_velocity_specific(const DHCellAccessor& dh_cell, arma::vec& solution,
552  double edge_scale, double edge_source_term)// override
553  {
554  const ElementAccessor<3> ele = dh_cell.elm();
555 
556  double storativity = ad_->storativity.value(ele.centre(), ele)
557  + ad_->extra_storativity.value(ele.centre(), ele);
558  double new_pressure, old_pressure, time_term = 0.0;
559 
560  for (unsigned int i=0; i<ele->n_sides(); i++) {
561 
562  if( ! ad_->use_steady_assembly_)
563  {
564  new_pressure = ad_->p_edge_solution.get(loc_schur_.row_dofs[i]);
565  old_pressure = ad_->p_edge_solution_previous_time.get(loc_schur_.row_dofs[i]);
566  time_term = edge_scale * storativity / ad_->time_step_ * (new_pressure - old_pressure);
567  }
568  solution[loc_side_dofs[i]] += edge_source_term - time_term;
569  }
570  }
571 
572  // assembly volume integrals
576 
577  NeighSideValues<dim<3?dim:2> ngh_values_;
578 
579  // Interpolation of velocity into barycenters
580  QGauss velocity_interpolation_quad_;
581  FEValues<3> velocity_interpolation_fv_;
582 
583  // data shared by assemblers of different dimension
584  AssemblyDataPtrLMH ad_;
585 
586  /** TODO: Investigate why the hell do we need this flag.
587  * If removed, it does not break any of the integration tests,
588  * however it must influence the Dirichlet rows in matrix.
589  */
590  std::vector<unsigned int> dirichlet_edge;
591 
592  LocalSystem loc_system_;
593  LocalSystem loc_schur_;
594  std::vector<unsigned int> loc_side_dofs;
595  std::vector<unsigned int> loc_edge_dofs;
596  unsigned int loc_ele_dof;
597 
598  // std::shared_ptr<MortarAssemblyBase> mortar_assembly;
599 
600  /// Index offset in the local system for the Schur complement.
601  unsigned int schur_offset_;
602 
603  /// Vector for reconstruted solution (velocity and pressure on element) from Schur complement.
604  arma::vec reconstructed_solution_;
605 
606  /// Flag for saving the local system.
607  /// Currently used only in case of seepage BC.
608  bool save_local_system_;
609 
610  /// Flag indicating whether the fluxes for seepage BC has been reconstructed already.
611  bool bc_fluxes_reconstruted;
612 };
613 
614 
615 #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:277
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:885
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:302
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:541
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
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:137
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.
Region region() const
Definition: accessors.hh:165
FE_RT0< dim > fe_rt_
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:79
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:66
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:223
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:85
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:296
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)