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