Flow123d  release_2.1.0-84-g6a13a75
darcy_flow_mh.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file darcy_flow_mh.cc
15  * @ingroup flow
16  * @brief Setup and solve linear system of mixed-hybrid discretization of the linear
17  * porous media flow with possible preferential flow in fractures and chanels.
18  */
19 
20 //#include <limits>
21 #include <vector>
22 //#include <iostream>
23 //#include <iterator>
24 //#include <algorithm>
25 #include <armadillo>
26 
27 #include "petscmat.h"
28 #include "petscviewer.h"
29 #include "petscerror.h"
30 #include "mpi.h"
31 
32 #include "system/system.hh"
33 #include "system/sys_profiler.hh"
34 #include "input/factory.hh"
35 
36 #include "mesh/mesh.h"
37 #include "mesh/intersection.hh"
38 #include "mesh/partitioning.hh"
39 #include "la/distribution.hh"
40 #include "la/linsys.hh"
41 #include "la/linsys_PETSC.hh"
42 #include "la/linsys_BDDC.hh"
43 #include "la/schur.hh"
44 //#include "la/sparse_graph.hh"
46 
47 #include "flow/darcy_flow_mh.hh"
48 
50 
51 /*
52 #include "fem/mapping_p1.hh"
53 #include "fem/fe_p.hh"
54 #include "fem/fe_values.hh"
55 #include "fem/fe_rt.hh"
56 #include "quadrature/quadrature_lib.hh"
57 */
58 
59 #include "tools/time_governor.hh"
61 #include "fields/field.hh"
62 #include "fields/field_values.hh"
64 
65 #include "coupling/balance.hh"
66 
67 #include "fields/vec_seq_double.hh"
68 #include "darcy_flow_assembly.hh"
69 
70 
71 FLOW123D_FORCE_LINK_IN_CHILD(darcy_flow_mh);
72 
73 
74 
75 
76 namespace it = Input::Type;
77 
79  return it::Selection("MH_MortarMethod")
80  .add_value(NoMortar, "None", "Mortar space: P0 on elements of lower dimension.")
81  .add_value(MortarP0, "P0", "Mortar space: P0 on elements of lower dimension.")
82  .add_value(MortarP1, "P1", "Mortar space: P1 on intersections, using non-conforming pressures.")
83  .close();
84 }
85 
86 
88  return it::Selection("Flow_Darcy_BC_Type")
89  .add_value(none, "none",
90  "Homogeneous Neumann boundary condition. Zero flux")
91  .add_value(dirichlet, "dirichlet",
92  "Dirichlet boundary condition. "
93  "Specify the pressure head through the ''bc_pressure'' field "
94  "or the piezometric head through the ''bc_piezo_head'' field.")
95  .add_value(total_flux, "total_flux", "Flux boundary condition (combines Neumann and Robin type). "
96  "Water inflow equal to (($q^N + \\sigma (h^R - h)$)). "
97  "Specify the water inflow by the 'bc_flux' field, the transition coefficient by 'bc_robin_sigma' "
98  "and the reference pressure head or pieozmetric head through ''bc_pressure'' or ''bc_piezo_head'' respectively.")
99  .add_value(seepage, "seepage",
100  "Seepage face boundary condition. Pressure and inflow bounded from above. Boundary with potential seepage flow "
101  "is described by the pair of inequalities: "
102  "(($h \\le h_d^D$)) and (($ q \\le q_d^N$)), where the equality holds in at least one of them. "
103  "Caution. Setting (($q_d^N$)) strictly negative "
104  "may lead to an ill posed problem since a positive outflow is enforced. "
105  "Parameters (($h_d^D$)) and (($q_d^N$)) are given by fields ``bc_pressure`` (or ``bc_piezo_head``) and ``bc_flux`` respectively."
106  )
107  .add_value(river, "river",
108  "River boundary condition. For the water level above the bedrock, (($H > H^S$)), the Robin boundary condition is used with the inflow given by: "
109  "(( $q^N + \\sigma(H^D - H)$)). For the water level under the bedrock, constant infiltration is used "
110  "(( $q^N + \\sigma(H^D - H^S)$)). Parameters: ``bc_pressure``, ``bc_switch_pressure``, "
111  " ``bc_sigma, ``bc_flux``."
112  )
113  .close();
114 }
115 
117 
118  const it::Record &field_descriptor =
119  it::Record("Flow_Darcy_MH_Data",FieldCommon::field_descriptor_record_description("Flow_Darcy_MH_Data") )
120  .copy_keys( DarcyMH::EqData().make_field_descriptor_type("Flow_Darcy_MH_Data_aux") )
121  .declare_key("bc_piezo_head", FieldAlgorithmBase< 3, FieldValue<3>::Scalar >::get_input_type_instance(),
122  "Boundary piezometric head for BC types: dirichlet, robin, and river." )
123  .declare_key("bc_switch_piezo_head", FieldAlgorithmBase< 3, FieldValue<3>::Scalar >::get_input_type_instance(),
124  "Boundary switch piezometric head for BC types: seepage, river." )
125  .declare_key("init_piezo_head", FieldAlgorithmBase< 3, FieldValue<3>::Scalar >::get_input_type_instance(),
126  "Initial condition for the pressure given as the piezometric head." )
127  .close();
128  return field_descriptor;
129 }
130 
132 
133  it::Record ns_rec = Input::Type::Record("NonlinearSolver", "Parameters to a non-linear solver.")
134  .declare_key("linear_solver", LinSys::get_input_type(), it::Default("{}"),
135  "Linear solver for MH problem.")
136  .declare_key("tolerance", it::Double(0.0), it::Default("1E-6"),
137  "Residual tolerance.")
138  .declare_key("min_it", it::Integer(0), it::Default("1"),
139  "Minimum number of iterations (linear solves) to use. This is usefull if the convergence criteria "
140  "does not characterize your goal well enough so it converges prematurely possibly without the single linear solve."
141  "If greater then 'max_it' the value is set to 'max_it'.")
142  .declare_key("max_it", it::Integer(0), it::Default("100"),
143  "Maximum number of iterations (linear solves) of the non-linear solver.")
144  .declare_key("converge_on_stagnation", it::Bool(), it::Default("false"),
145  "If a stagnation of the nonlinear solver is detected the solver stops. "
146  "A divergence is reported by default forcing the end of the simulation. Setting this flag to 'true', the solver"
147  "ends with convergence success on stagnation, but report warning about it.")
148  .close();
149 
150  return it::Record("Flow_Darcy_MH", "Mixed-Hybrid solver for STEADY saturated Darcy flow.")
152  .declare_key("gravity", it::Array(it::Double(), 3,3), it::Default("[ 0, 0, -1]"),
153  "Vector of the gravity force. Dimensionless.")
155  "Input data for Darcy flow model.")
156  .declare_key("nonlinear_solver", ns_rec, it::Default::obligatory(),
157  "Non-linear solver for MH problem.")
159  "Parameters of output stream.")
160 
161  .declare_key("output", DarcyFlowMHOutput::get_input_type(), IT::Default("{ \"fields\": [ \"pressure_p0\", \"velocity_p0\" ] }"),
162  "Parameters of output from MH module.")
164  "Parameters of output form MH module.")
165  .declare_key("balance", Balance::get_input_type(), it::Default("{}"),
166  "Settings for computing mass balance.")
168  "Time governor setting for the unsteady Darcy flow model.")
169  .declare_key("n_schurs", it::Integer(0,2), it::Default("2"),
170  "Number of Schur complements to perform when solving MH system.")
171  .declare_key("mortar_method", get_mh_mortar_selection(), it::Default("\"None\""),
172  "Method for coupling Darcy flow between dimensions." )
173  .close();
174 }
175 
176 
177 const int DarcyMH::registrar =
178  Input::register_class< DarcyMH, Mesh &, const Input::Record >("Flow_Darcy_MH") +
180 
181 
182 
184 {
185  ADD_FIELD(anisotropy, "Anisotropy of the conductivity tensor.", "1.0" );
186  anisotropy.units( UnitSI::dimensionless() );
187 
188  ADD_FIELD(cross_section, "Complement dimension parameter (cross section for 1D, thickness for 2D).", "1.0" );
189  cross_section.units( UnitSI().m(3).md() );
190 
191  ADD_FIELD(conductivity, "Isotropic conductivity scalar.", "1.0" );
192  conductivity.units( UnitSI().m().s(-1) ).set_limits(0.0);
193 
194  ADD_FIELD(sigma, "Transition coefficient between dimensions.", "1.0");
195  sigma.units( UnitSI::dimensionless() );
196 
197  ADD_FIELD(water_source_density, "Water source density.", "0.0");
198  water_source_density.units( UnitSI().s(-1) );
199 
200  ADD_FIELD(bc_type,"Boundary condition type, possible values:", "\"none\"" );
201  // TODO: temporary solution, we should try to get rid of pointer to the selection after having generic types
202  bc_type.input_selection( get_bc_type_selection() );
203  bc_type.units( UnitSI::dimensionless() );
204 
205  ADD_FIELD(bc_pressure,"Prescribed pressure value on the boundary. Used for all values of 'bc_type' save the bc_type='none'."
206  "See documentation of 'bc_type' for exact meaning of 'bc_pressure' in individual boundary condition types.", "0.0");
207  bc_pressure.disable_where(bc_type, {none} );
208  bc_pressure.units( UnitSI().m() );
209 
210  ADD_FIELD(bc_flux,"Incoming water boundary flux. Used for bc_types : 'none', 'total_flux', 'seepage', 'river'.", "0.0");
211  bc_flux.disable_where(bc_type, {none, dirichlet} );
212  bc_flux.units( UnitSI().m(4).s(-1).md() );
213 
214  ADD_FIELD(bc_robin_sigma,"Conductivity coefficient in the 'total_flux' or the 'river' boundary condition type.", "0.0");
215  bc_robin_sigma.disable_where(bc_type, {none, dirichlet, seepage} );
216  bc_robin_sigma.units( UnitSI().m(3).s(-1).md() );
217 
218  ADD_FIELD(bc_switch_pressure,
219  "Critical switch pressure for 'seepage' and 'river' boundary conditions.", "0.0");
220  bc_switch_pressure.disable_where(bc_type, {none, dirichlet, total_flux} );
221  bc_switch_pressure.units( UnitSI().m() );
222 
223  //these are for unsteady
224  ADD_FIELD(init_pressure, "Initial condition as pressure", "0.0" );
225  init_pressure.units( UnitSI().m() );
226 
227  ADD_FIELD(storativity,"Storativity.", "0.0" );
228  storativity.units( UnitSI().m(-1) );
229 
230  //time_term_fields = this->subset({"storativity"});
231  //main_matrix_fields = this->subset({"anisotropy", "conductivity", "cross_section", "sigma", "bc_type", "bc_robin_sigma"});
232  //rhs_fields = this->subset({"water_source_density", "bc_pressure", "bc_flux"});
233 
234 }
235 
236 
237 
238 
239 
240 
241 //=============================================================================
242 // CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
243 // - do it in parallel:
244 // - initial distribution of elements, edges
245 //
246 /*! @brief CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
247  *
248  * Parameters {Solver,NSchurs} number of performed Schur
249  * complements (0,1,2) for water flow MH-system
250  *
251  */
252 //=============================================================================
253 DarcyMH::DarcyMH(Mesh &mesh_in, const Input::Record in_rec)
254 : DarcyFlowInterface(mesh_in, in_rec),
255  solution(nullptr),
256  schur0(nullptr),
257  data_changed_(false),
258  output_object(nullptr)
259 {
260 
261  is_linear_=true;
262 
263  START_TIMER("Darcy constructor");
264  {
265  auto time_record = input_record_.val<Input::Record>("time");
266  //if ( in_rec.opt_val("time", time_record) )
267  time_ = new TimeGovernor(time_record);
268  //else
269  // time_ = new TimeGovernor();
270  }
271 
272  data_ = make_shared<EqData>();
274 
275 
276  size = mesh_->n_elements() + mesh_->n_sides() + mesh_->n_edges();
277  n_schur_compls = in_rec.val<int>("n_schurs");
278  mortar_method_= in_rec.val<MortarMethod>("mortar_method");
279  if (mortar_method_ != NoMortar) {
281  }
282 
283 
284 
285  //side_ds->view( std::cout );
286  //el_ds->view( std::cout );
287  //edge_ds->view( std::cout );
288  //rows_ds->view( std::cout );
289 
290 }
291 
292 
293 
295 //connecting data fields with mesh
296 {
297 
298  START_TIMER("data init");
299  data_->mesh = mesh_;
300  data_->mh_dh = &mh_dh;
301  data_->set_mesh(*mesh_);
302 
303  auto gravity_array = input_record_.val<Input::Array>("gravity");
304  std::vector<double> gvec;
305  gravity_array.copy_to(gvec);
306  gvec.push_back(0.0); // zero pressure shift
307  data_->gravity_ = arma::vec(gvec);
308  data_->gravity_vec_ = data_->gravity_.subvec(0,2);
309 
310  data_->bc_pressure.add_factory(
311  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
312  (data_->gravity_, "bc_piezo_head") );
313  data_->bc_switch_pressure.add_factory(
314  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
315  (data_->gravity_, "bc_switch_piezo_head") );
316  data_->init_pressure.add_factory(
317  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
318  (data_->gravity_, "init_piezo_head") );
319 
320 
321  data_->set_input_list( this->input_record_.val<Input::Array>("input_fields") );
322  // Check that the time step was set for the transient simulation.
323  if (! zero_time_term(true) && time_->is_default() ) {
324  //THROW(ExcAssertMsg());
325  //THROW(ExcMissingTimeGovernor() << input_record_.ei_address());
326  MessageOut() << "Missing the key 'time', obligatory for the transient problems." << endl;
327  ASSERT(false);
328  }
329 
330  data_->mark_input_times(*time_);
331 }
332 
333 
334 
336 
337  init_eq_data();
339 
340  mh_dh.reinit(mesh_);
341  // Initialize bc_switch_dirichlet to size of global boundary.
343 
344 
347  .val<Input::Record>("nonlinear_solver")
348  .val<Input::AbstractRecord>("linear_solver");
349 
350  // auxiliary set_time call since allocation assembly evaluates fields as well
352  data_->system_.local_matrix = std::make_shared<arma::mat>();
354 
355 
356 
357  // allocate time term vectors
358  VecDuplicate(schur0->get_solution(), &previous_solution);
359  VecCreateMPI(PETSC_COMM_WORLD, mh_dh.rows_ds->lsize(),PETSC_DETERMINE,&(steady_diagonal));
360  VecDuplicate(steady_diagonal,& new_diagonal);
361  VecZeroEntries(new_diagonal);
362  VecDuplicate(steady_diagonal, &steady_rhs);
363 
364 
365  // initialization of balance object
366  balance_ = std::make_shared<Balance>("water", mesh_);
367  balance_->init_from_input(input_record_.val<Input::Record>("balance"), time());
368  if (balance_)
369  {
370  data_-> water_balance_idx_ = water_balance_idx_ = balance_->add_quantity("water_volume");
371  balance_->allocate(mh_dh.rows_ds->lsize(), 1);
372  balance_->units(UnitSI().m(3));
373  }
374 
375 
376  data_->system_.balance = balance_;
377  data_->system_.lin_sys = schur0;
378 
379 
381 }
382 
384 {}
385 
387 {
388 
389  /* TODO:
390  * - Allow solution reconstruction (pressure and velocity) from initial condition on user request.
391  * - Steady solution as an intitial condition may be forced by setting inti_time =-1, and set data for the steady solver in that time.
392  * Solver should be able to switch from and to steady case depending on the zero time term.
393  */
394 
396 
397  // zero_time_term means steady case
398  bool zero_time_term_from_right = zero_time_term();
399 
400 
401  if (zero_time_term_from_right) {
402  // steady case
403  VecZeroEntries(schur0->get_solution());
404  //read_initial_condition(); // Possible solution guess for steady case.
405  use_steady_assembly_ = true;
406  solve_nonlinear(); // with right limit data
407  } else {
408  VecZeroEntries(schur0->get_solution());
409  VecZeroEntries(previous_solution);
411  assembly_linear_system(); // in particular due to balance
412  // TODO: reconstruction of solution in zero time.
413  }
414  //solution_output(T,right_limit); // data for time T in any case
415  output_data();
416 }
417 
418 //=============================================================================
419 // COMPOSE and SOLVE WATER MH System possibly through Schur complements
420 //=============================================================================
422 {
423  START_TIMER("Solving MH system");
424 
425  time_->next_time();
426 
427  time_->view("DARCY"); //time governor information output
429  bool zero_time_term_from_left=zero_time_term();
430 
431  bool jump_time = data_->storativity.is_jump_time();
432  if (! zero_time_term_from_left) {
433  // time term not treated as zero
434  // Unsteady solution up to the T.
435 
436  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
437  use_steady_assembly_ = false;
438  prepare_new_time_step(); //SWAP
439 
440  solve_nonlinear(); // with left limit data
441  if (jump_time) {
442  WarningOut() << "Output of solution discontinuous in time not supported yet.\n";
443  //solution_output(T, left_limit); // output use time T- delta*dt
444  //output_data();
445  }
446  }
447 
448  if (time_->is_end()) {
449  // output for unsteady case, end_time should not be the jump time
450  // but rether check that
451  if (! zero_time_term_from_left && ! jump_time) output_data();
452  return;
453  }
454 
456  bool zero_time_term_from_right=zero_time_term();
457  if (zero_time_term_from_right) {
458  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
459  use_steady_assembly_ = true;
460  solve_nonlinear(); // with right limit data
461 
462  } else if (! zero_time_term_from_left && jump_time) {
463  WarningOut() << "Discontinuous time term not supported yet.\n";
464  //solution_transfer(); // internally call set_time(T, left) and set_time(T,right) again
465  //solve_nonlinear(); // with right limit data
466  }
467  //solution_output(T,right_limit); // data for time T in any case
468  output_data();
469 
470 }
471 
472 bool DarcyMH::zero_time_term(bool time_global) {
473  if (time_global) {
474  return (data_->storativity.input_list_size() == 0);
475  } else {
476  return data_->storativity.field_result(mesh_->region_db().get_region_set("BULK")) == result_zeros;
477  }
478 }
479 
480 
482 {
483 
485  double residual_norm = schur0->compute_residual();
486  unsigned int l_it=0;
488  MessageOut().fmt("[nonlin solver] norm of initial residual: {}\n", residual_norm);
489 
490  // Reduce is_linear flag.
491  int is_linear_common;
492  MPI_Allreduce(&is_linear_, &is_linear_common,1, MPI_INT ,MPI_MIN,PETSC_COMM_WORLD);
493 
494  Input::Record nl_solver_rec = input_record_.val<Input::Record>("nonlinear_solver");
495  this->tolerance_ = nl_solver_rec.val<double>("tolerance");
496  this->max_n_it_ = nl_solver_rec.val<unsigned int>("max_it");
497  this->min_n_it_ = nl_solver_rec.val<unsigned int>("min_it");
498  if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->max_n_it_;
499 
500  if (! is_linear_common) {
501  // set tolerances of the linear solver unless they are set by user.
502  schur0->set_tolerances(0.1*this->tolerance_, 0.01*this->tolerance_, 100);
503  }
504  vector<double> convergence_history;
505 
506  Vec save_solution;
507  VecDuplicate(schur0->get_solution(), &save_solution);
508  while (nonlinear_iteration_ < this->min_n_it_ ||
509  (residual_norm > this->tolerance_ && nonlinear_iteration_ < this->max_n_it_ )) {
510  OLD_ASSERT_EQUAL( convergence_history.size(), nonlinear_iteration_ );
511  convergence_history.push_back(residual_norm);
512 
513  // stagnation test
514  if (convergence_history.size() >= 5 &&
515  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
516  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
517  // stagnation
518  if (input_record_.val<Input::Record>("nonlinear_solver").val<bool>("converge_on_stagnation")) {
519  WarningOut().fmt("Accept solution on stagnation. Its: {} Residual: {}\n", nonlinear_iteration_, residual_norm);
520  break;
521  } else {
522  THROW(ExcSolverDiverge() << EI_Reason("Stagnation."));
523  }
524  }
525 
526  if (! is_linear_common)
527  VecCopy( schur0->get_solution(), save_solution);
528  int convergedReason = schur0->solve();
530 
531  // hack to make BDDC work with empty compute_residual
532  if (is_linear_common) break;
533  data_changed_=true; // force reassembly for non-linear case
534 
535  double alpha = 1; // how much of new solution
536  VecAXPBY(schur0->get_solution(), (1-alpha), alpha, save_solution);
537 
538  /*
539  double * sol;
540  unsigned int sol_size;
541  get_solution_vector(sol, sol_size);
542  if (mh_dh.el_ds->myp() == 0)
543  VecView(sol_vec, PETSC_VIEWER_STDOUT_SELF);
544  */
545 
546  //LogOut().fmt("Linear solver ended with reason: {} \n", convergedReason );
547  //OLD_ASSERT( convergedReason >= 0, "Linear solver failed to converge. Convergence reason %d \n", convergedReason );
549 
550  residual_norm = schur0->compute_residual();
551  MessageOut().fmt("[nonlinear solver] it: {} lin. it:{} (reason: {}) residual: {}\n",
552  nonlinear_iteration_, l_it, convergedReason, residual_norm);
553  }
554  this -> postprocess();
555 
556  // adapt timestep
557  if (! this->zero_time_term()) {
558  double mult = 1.0;
559  if (nonlinear_iteration_ < 3) mult = 1.6;
560  if (nonlinear_iteration_ > 7) mult = 0.7;
561  int result = time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
562  //DebugOut().fmt("time adaptivity, res: {} it: {} m: {} dt: {} edt: {}\n", result, nonlinear_iteration_, mult, time_->dt(), time_->estimate_dt());
563  }
564 
566 
567 }
568 
569 
571 {
572  //VecSwap(previous_solution, schur0->get_solution());
573 }
574 
576 {
577  START_TIMER("postprocess");
578  //ElementFullIter ele = ELEMENT_FULL_ITER(mesh_, NULL);
579 
580  // modify side fluxes in parallel
581  // for every local edge take time term on digonal and add it to the corresponding flux
582  /*
583  for (unsigned int i_loc = 0; i_loc < el_ds->lsize(); i_loc++) {
584  ele = mesh_->element(el_4_loc[i_loc]);
585  FOR_ELEMENT_SIDES(ele,i) {
586  side_rows[i] = side_row_4_id[ mh_dh.side_dof( ele_ac.side(i) ) ];
587  values[i] = -1.0 * ele_ac.measure() *
588  data.cross_section.value(ele_ac.centre(), ele_ac.element_accessor()) *
589  data.water_source_density.value(ele_ac.centre(), ele_ac.element_accessor()) /
590  ele_ac.n_sides();
591  }
592  VecSetValues(schur0->get_solution(), ele_ac.n_sides(), side_rows, values, ADD_VALUES);
593  }
594  VecAssemblyBegin(schur0->get_solution());
595  VecAssemblyEnd(schur0->get_solution());
596  */
597 }
598 
599 
601  START_TIMER("Darcy output data");
602  //time_->view("DARCY"); //time governor information output
603  this->output_object->output();
604 
605 
606  if (balance_ != nullptr)
607  {
608  START_TIMER("Darcy balance output");
609  if (balance_->cumulative() && time_->tlevel() > 0)
610  {
611  balance_->calculate_cumulative_sources(water_balance_idx_, schur0->get_solution(), time_->dt());
612  balance_->calculate_cumulative_fluxes(water_balance_idx_, schur0->get_solution(), time_->dt());
613  }
614 
615  if ( balance_->is_current( time().step()) )
616  {
617  balance_->calculate_mass(water_balance_idx_, schur0->get_solution());
618  balance_->calculate_source(water_balance_idx_, schur0->get_solution());
619  balance_->calculate_flux(water_balance_idx_, schur0->get_solution());
620  balance_->output(time_->t());
621  }
622  }
623 }
624 
625 
627 {
628  return schur0->get_solution_precision();
629 }
630 
631 
632 
633 void DarcyMH::get_solution_vector(double * &vec, unsigned int &vec_size)
634 {
635  // TODO: make class for vectors (wrapper for PETSC or other) derived from LazyDependency
636  // and use its mechanism to manage dependency between vectors
638 
639  // scatter solution to all procs
640  VecScatterBegin(par_to_all, schur0->get_solution(), sol_vec, INSERT_VALUES, SCATTER_FORWARD);
641  VecScatterEnd( par_to_all, schur0->get_solution(), sol_vec, INSERT_VALUES, SCATTER_FORWARD);
643  }
644 
645  vec = solution;
646  vec_size = this->size;
647  OLD_ASSERT(vec != NULL, "Requested solution is not allocated!\n");
648 }
649 
651 {
652  vec=schur0->get_solution();
653  OLD_ASSERT(vec != NULL, "Requested solution is not allocated!\n");
654 }
655 
656 
657 /*
658 void DarcyMH::local_assembly_specific() {
659 
660 }
661 */
662 // ===========================================================================================
663 //
664 // MATRIX ASSEMBLY - we use abstract assembly routine, where LS Mat/Vec SetValues
665 // are in fact pointers to allocating or filling functions - this is governed by Linsystem roitunes
666 //
667 // =======================================================================================
668 
670 {
671  START_TIMER("DarcyFlowMH_Steady::assembly_steady_mh_matrix");
672 
673  LinSys *ls = schur0;
674 
675  class Boundary *bcd;
676  class Neighbour *ngh;
677 
678  bool fill_matrix = assembler.size() > 0;
679  int side_row, edge_row, loc_b = 0;
680  int tmp_rows[100];
681  double local_vb[4]; // 2x2 matrix
682  int side_rows[4], edge_rows[4];
683 
684  // to make space for second schur complement, max. 10 neighbour edges of one el.
685  double zeros[1000];
686  for(int i=0; i<1000; i++) zeros[i]=0.0;
687 
688  double minus_ones[4] = { -1.0, -1.0, -1.0, -1.0 };
689  double * loc_side_rhs = data_->system_.loc_side_rhs;
690 
691  arma::mat &local_matrix = *(data_->system_.local_matrix);
692 
693 
694  if (balance_ != nullptr && fill_matrix)
695  balance_->start_flux_assembly(water_balance_idx_);
696 
697  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
698  auto ele_ac = mh_dh.accessor(i_loc);
699  unsigned int nsides = ele_ac.n_sides();
700  data_->system_.dirichlet_edge.resize(nsides);
701 
702 
703 
704  for (unsigned int i = 0; i < nsides; i++) {
705 
706 /* if (! side_ds->is_local(idx_side)) {
707  cout << el_ds->myp() << " : iside: " << ele.index() << " [" << el_ds->begin() << ", " << el_ds->end() << "]" << endl;
708  cout << el_ds->myp() << " : iside: " << idx_side << " [" << side_ds->begin() << ", " << side_ds->end() << "]" << endl;
709 
710  }*/
711 
712  side_rows[i] = side_row = ele_ac.side_row(i);
713  edge_rows[i] = edge_row = ele_ac.edge_row(i);
714  bcd=ele_ac.side(i)->cond();
715 
716  // gravity term on RHS
717  //
718  loc_side_rhs[i] = 0;
719 
720  // set block C and C': side-edge, edge-side
721  double c_val = 1.0;
722  data_->system_.dirichlet_edge[i] = 0;
723 
724  if (bcd) {
725  ElementAccessor<3> b_ele = bcd->element_accessor();
726  EqData::BC_Type type = (EqData::BC_Type)data_->bc_type.value(b_ele.centre(), b_ele);
727 
728  double cross_section = data_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor());
729 
730  if ( type == EqData::none) {
731  // homogeneous neumann
732  } else if ( type == EqData::dirichlet ) {
733  double bc_pressure = data_->bc_pressure.value(b_ele.centre(), b_ele);
734  c_val = 0.0;
735  loc_side_rhs[i] -= bc_pressure;
736  ls->rhs_set_value(edge_row, -bc_pressure);
737  ls->mat_set_value(edge_row, edge_row, -1.0);
738  data_->system_.dirichlet_edge[i] = 1;
739 
740  } else if ( type == EqData::total_flux) {
741 
742  // internally we work with outward flux
743  double bc_flux = -data_->bc_flux.value(b_ele.centre(), b_ele);
744  double bc_pressure = data_->bc_pressure.value(b_ele.centre(), b_ele);
745  double bc_sigma = data_->bc_robin_sigma.value(b_ele.centre(), b_ele);
746  //DebugOut().fmt("erow: {} flux: {} mesure: {} cs: {}", edge_row, bc_flux, bcd->element()->measure(), cross_section);
747  ls->mat_set_value(edge_row, edge_row, -bcd->element()->measure() * bc_sigma * cross_section );
748  ls->rhs_set_value(edge_row, (bc_flux - bc_sigma * bc_pressure) * bcd->element()->measure() * cross_section);
749 
750  } else if (type==EqData::seepage) {
751  is_linear_=false;
752  //unsigned int loc_edge_idx = edge_row - rows_ds->begin() - side_ds->lsize() - el_ds->lsize();
753  unsigned int loc_edge_idx = bcd->bc_ele_idx_;
754  char & switch_dirichlet = bc_switch_dirichlet[loc_edge_idx];
755  double bc_pressure = data_->bc_switch_pressure.value(b_ele.centre(), b_ele);
756  double bc_flux = -data_->bc_flux.value(b_ele.centre(), b_ele);
757  double side_flux=bc_flux * bcd->element()->measure() * cross_section;
758 
759  // ** Update BC type. **
760  if (switch_dirichlet) {
761  // check and possibly switch to flux BC
762  // The switch raise error on the corresponding edge row.
763  // Magnitude of the error is abs(solution_flux - side_flux).
764  ASSERT_DBG(mh_dh.rows_ds->is_local(side_row))(side_row);
765  unsigned int loc_side_row = ele_ac.side_local_row(i);
766  double & solution_flux = ls->get_solution_array()[loc_side_row];
767 
768  if ( solution_flux < side_flux) {
769  //DebugOut().fmt("x: {}, to neum, p: {} f: {} -> f: {}\n", b_ele.centre()[0], bc_pressure, solution_flux, side_flux);
770  solution_flux = side_flux;
771  switch_dirichlet=0;
772 
773  }
774  } else {
775  // check and possibly switch to pressure BC
776  // TODO: What is the appropriate DOF in not local?
777  // The switch raise error on the corresponding side row.
778  // Magnitude of the error is abs(solution_head - bc_pressure)
779  // Since usually K is very large, this error would be much
780  // higher then error caused by the inverse switch, this
781  // cause that a solution with the flux violating the
782  // flux inequality leading may be accepted, while the error
783  // in pressure inequality is always satisfied.
784  ASSERT_DBG(mh_dh.rows_ds->is_local(edge_row))(edge_row);
785  unsigned int loc_edge_row = ele_ac.edge_local_row(i);
786  double & solution_head = ls->get_solution_array()[loc_edge_row];
787 
788  if ( solution_head > bc_pressure) {
789  //DebugOut().fmt("x: {}, to dirich, p: {} -> p: {} f: {}\n",b_ele.centre()[0], solution_head, bc_pressure, bc_flux);
790  solution_head = bc_pressure;
791  switch_dirichlet=1;
792  }
793  }
794 
795  // ** Apply BCUpdate BC type. **
796  // Force Dirichlet type during the first iteration of the unsteady case.
797  if (switch_dirichlet || (use_steady_assembly_ && nonlinear_iteration_ == 0) ) {
798  //DebugOut().fmt("x: {}, dirich, bcp: {}\n", b_ele.centre()[0], bc_pressure);
799  c_val = 0.0;
800  loc_side_rhs[i] -= bc_pressure;
801  ls->rhs_set_value(edge_row, -bc_pressure);
802  ls->mat_set_value(edge_row, edge_row, -1.0);
803  data_->system_.dirichlet_edge[i] = 1;
804  } else {
805  //DebugOut()("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], side_flux, bc_flux);
806  ls->rhs_set_value(edge_row, side_flux);
807  }
808 
809  } else if (type==EqData::river) {
810  is_linear_=false;
811  //unsigned int loc_edge_idx = edge_row - rows_ds->begin() - side_ds->lsize() - el_ds->lsize();
812  //unsigned int loc_edge_idx = bcd->bc_ele_idx_;
813  //char & switch_dirichlet = bc_switch_dirichlet[loc_edge_idx];
814 
815  double bc_pressure = data_->bc_pressure.value(b_ele.centre(), b_ele);
816  double bc_switch_pressure = data_->bc_switch_pressure.value(b_ele.centre(), b_ele);
817  double bc_flux = -data_->bc_flux.value(b_ele.centre(), b_ele);
818  double bc_sigma = data_->bc_robin_sigma.value(b_ele.centre(), b_ele);
819  ASSERT_DBG(mh_dh.rows_ds->is_local(edge_row))(edge_row);
820  unsigned int loc_edge_row = ele_ac.edge_local_row(i);
821  double & solution_head = ls->get_solution_array()[loc_edge_row];
822 
823 
824  // Force Robin type during the first iteration of the unsteady case.
825  if (solution_head > bc_switch_pressure || (use_steady_assembly_ && nonlinear_iteration_ ==0)) {
826  // Robin BC
827  //DebugOut().fmt("x: {}, robin, bcp: {}\n", b_ele.centre()[0], bc_pressure);
828  ls->rhs_set_value(edge_row, bcd->element()->measure() * cross_section * (bc_flux - bc_sigma * bc_pressure) );
829  ls->mat_set_value(edge_row, edge_row, -bcd->element()->measure() * bc_sigma * cross_section );
830  } else {
831  // Neumann BC
832  //DebugOut().fmt("x: {}, neuman, q: {} bcq: {}\n", b_ele.centre()[0], bc_switch_pressure, bc_pressure);
833  double bc_total_flux = bc_flux + bc_sigma*(bc_switch_pressure - bc_pressure);
834  ls->rhs_set_value(edge_row, bc_total_flux * bcd->element()->measure() * cross_section);
835  }
836  } else {
837  xprintf(UsrErr, "BC type not supported.\n");
838  }
839 
840  if (balance_ != nullptr && fill_matrix)
841  {
842  /*
843  DebugOut()("add_flux: {} {} {} {}\n",
844  mh_dh.el_ds->myp(),
845  ele_ac.ele_global_idx(),
846  loc_b,
847  side_row);*/
848  balance_->add_flux_matrix_values(water_balance_idx_, loc_b, {side_row}, {1});
849  }
850  ++loc_b;
851  }
852  ls->mat_set_value(side_row, edge_row, c_val);
853  ls->mat_set_value(edge_row, side_row, c_val);
854 
855  }
856 
857  if (fill_matrix) {
858  assembler[ele_ac.dim()-1]->assembly_local_matrix(ele_ac);
859 
860  // assemble matrix for weights in BDDCML
861  // approximation to diagonal of
862  // S = -C - B*inv(A)*B'
863  // as
864  // diag(S) ~ - diag(C) - 1./diag(A)
865  // the weights form a partition of unity to average a discontinuous solution from neighbouring subdomains
866  // to a continuous one
867  // it is important to scale the effect - if conductivity is low for one subdomain and high for the other,
868  // trust more the one with low conductivity - it will be closer to the truth than an arithmetic average
869  if ( typeid(*ls) == typeid(LinSys_BDDC) ) {
870  for(unsigned int i=0; i < nsides; i++) {
871  double val_side = local_matrix(i,i);
872  double val_edge = -1./local_matrix(i,i);
873 
874  static_cast<LinSys_BDDC*>(ls)->diagonal_weights_set_value( side_rows[i], val_side );
875  static_cast<LinSys_BDDC*>(ls)->diagonal_weights_set_value( edge_rows[i], val_edge );
876  }
877  }
878  }
879 
880  ls->rhs_set_values(nsides, side_rows, loc_side_rhs);
881 
882 
883  // set block A: side-side on one element - block diagonal matrix
884  ls->mat_set_values(nsides, side_rows, nsides, side_rows, local_matrix.memptr());
885  // set block B, B': element-side, side-element
886  int ele_row = ele_ac.ele_row();
887  ls->mat_set_values(1, &ele_row, nsides, side_rows, minus_ones);
888  ls->mat_set_values(nsides, side_rows, 1, &ele_row, minus_ones);
889 
890 
891  // D block: diagonal: element-element
892 
893  ls->mat_set_value(ele_row, ele_row, 0.0); // maybe this should be in virtual block for schur preallocation
894 
895  if ( typeid(*ls) == typeid(LinSys_BDDC) ) {
896  double val_ele = 1.;
897  static_cast<LinSys_BDDC*>(ls)->diagonal_weights_set_value( ele_row, val_ele );
898  }
899 
900  // D, E',E block: compatible connections: element-edge
901 
902  for (unsigned int i = 0; i < ele_ac.full_iter()->n_neighs_vb; i++) {
903  // every compatible connection adds a 2x2 matrix involving
904  // current element pressure and a connected edge pressure
905  ngh= ele_ac.full_iter()->neigh_vb[i];
906  tmp_rows[0]=ele_row;
907  tmp_rows[1]=mh_dh.row_4_edge[ ngh->edge_idx() ];
908 
909  if (fill_matrix)
910  assembler[ngh->side()->dim()]->assembly_local_vb(local_vb, ele_ac.full_iter(), ngh);
911 
912  ls->mat_set_values(2, tmp_rows, 2, tmp_rows, local_vb);
913 
914  // update matrix for weights in BDDCML
915  if ( typeid(*ls) == typeid(LinSys_BDDC) ) {
916  int ind = tmp_rows[1];
917  // there is -value on diagonal in block C!
918  double new_val = local_vb[0];
919  static_cast<LinSys_BDDC*>(ls)->diagonal_weights_set_value( ind, new_val );
920  }
921 
922  if (n_schur_compls == 2) {
923  // for 2. Schur: N dim edge is conected with N dim element =>
924  // there are nz between N dim edge and N-1 dim edges of the element
925 
926  ls->mat_set_values(nsides, edge_rows, 1, tmp_rows+1, zeros);
927  ls->mat_set_values(1, tmp_rows+1, nsides, edge_rows, zeros);
928 
929  // save all global edge indices to higher positions
930  tmp_rows[2+i] = tmp_rows[1];
931  }
932  }
933 
934 
935  // add virtual values for schur complement allocation
936  uint n_neigh;
937  switch (n_schur_compls) {
938  case 2:
939  n_neigh = ele_ac.full_iter()->n_neighs_vb;
940  // Connections between edges of N+1 dim. elements neighboring with actual N dim element 'ele'
941  OLD_ASSERT(n_neigh*n_neigh<1000, "Too many values in E block.");
942  ls->mat_set_values(ele_ac.full_iter()->n_neighs_vb, tmp_rows+2,
943  ele_ac.full_iter()->n_neighs_vb, tmp_rows+2, zeros);
944 
945  case 1: // included also for case 2
946  // -(C')*(A-)*B block and its transpose conect edge with its elements
947  ls->mat_set_values(1, &ele_row, ele_ac.n_sides(), edge_rows, zeros);
948  ls->mat_set_values(ele_ac.n_sides(), edge_rows, 1, &ele_row, zeros);
949  // -(C')*(A-)*C block conect all edges of every element
950  ls->mat_set_values(ele_ac.n_sides(), edge_rows, ele_ac.n_sides(), edge_rows, zeros);
951  }
952  }
953 
954  if (balance_ != nullptr && fill_matrix)
955  balance_->finish_flux_assembly(water_balance_idx_);
956 
957 
958 
959  if (mortar_method_ == MortarP0) {
960  P0_CouplingAssembler(*this).assembly(*ls);
961  } else if (mortar_method_ == MortarP1) {
962  P1_CouplingAssembler(*this).assembly(*ls);
963  }
964 }
965 
966 
968 {
969  START_TIMER("assembly source term");
970  if (balance_ != nullptr)
971  balance_->start_source_assembly(water_balance_idx_);
972 
973  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
974  auto ele_ac = mh_dh.accessor(i_loc);
975 
976  double cs = data_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor());
977 
978  // set sources
979  double source = ele_ac.measure() * cs *
980  data_->water_source_density.value(ele_ac.centre(), ele_ac.element_accessor());
981  schur0->rhs_set_value(ele_ac.ele_row(), -1.0 * source );
982 
983  if (balance_ != nullptr)
984  balance_->add_source_vec_values(water_balance_idx_, ele_ac.region().bulk_idx(), {(int) ele_ac.ele_row()}, {source});
985  }
986 
987  if (balance_ != nullptr)
988  balance_->finish_source_assembly(water_balance_idx_);
989 }
990 
991 
993  vector<int> &dofs, unsigned int &ele_type, double &delta, arma::vec &dirichlet) {
994 
995  const Element *ele;
996 
997  if (i_ele == (int)(ml_it_->size()) ) { // master element .. 1D
998  ele_type = 0;
999  delta = -delta_0;
1000  ele=master_;
1001  } else {
1002  ele_type = 1;
1003  const Intersection &isect=intersections_[ (*ml_it_)[i_ele] ];
1004  delta = isect.intersection_true_size();
1005  ele = isect.slave_iter();
1006  }
1007 
1008  dofs.resize(ele->n_sides());
1009  dirichlet.resize(ele->n_sides());
1010  dirichlet.zeros();
1011 
1012  for(unsigned int i_side=0; i_side < ele->n_sides(); i_side++ ) {
1013  dofs[i_side]=darcy_.mh_dh.row_4_edge[ele->side(i_side)->edge_idx()];
1014  Boundary * bcd = ele->side(i_side)->cond();
1015  if (bcd) {
1016  ElementAccessor<3> b_ele = bcd->element_accessor();
1017  auto type = (DarcyMH::EqData::BC_Type)darcy_.data_->bc_type.value(b_ele.centre(), b_ele);
1018  //DebugOut().fmt("bcd id: {} sidx: {} type: {}\n", ele->id(), i_side, type);
1019  if (type == DarcyMH::EqData::dirichlet) {
1020  //DebugOut().fmt("Dirichlet: {}\n", ele->index());
1021  dofs[i_side] = -dofs[i_side];
1022  double bc_pressure = darcy_.data_->bc_pressure.value(b_ele.centre(), b_ele);
1023  dirichlet[i_side] = bc_pressure;
1024  }
1025  }
1026  }
1027 
1028 }
1029 
1030 /**
1031  * Works well but there is large error next to the boundary.
1032  */
1034  double delta_i, delta_j;
1035  arma::mat product;
1036  arma::vec dirichlet_i, dirichlet_j;
1037  unsigned int ele_type_i, ele_type_j; // element type 0-master, 1-slave for row and col
1038 
1039  unsigned int i,j;
1040  vector<int> dofs_i,dofs_j;
1041 
1042  for(ml_it_ = master_list_.begin(); ml_it_ != master_list_.end(); ++ml_it_) {
1043 
1044  if (ml_it_->size() == 0) continue; // skip empty masters
1045 
1046 
1047  // on the intersection element we consider
1048  // * intersection dofs for master and slave
1049  // those are dofs of the space into which we interpolate
1050  // base functions from individual master and slave elements
1051  // For the master dofs both are usualy eqivalent.
1052  // * original dofs - for master same as intersection dofs, for slave
1053  // all dofs of slave elements
1054 
1055  // form list of intersection dofs, in this case pressures in barycenters
1056  // but we do not use those form MH system in order to allow second schur somplement. We rather map these
1057  // dofs to pressure traces, i.e. we use averages of traces as barycentric values.
1058 
1059 
1060  master_ = intersections_[ml_it_->front()].master_iter();
1061  delta_0 = master_->measure();
1062 
1063  double master_sigma=darcy_.data_->sigma.value( master_->centre(), master_->element_accessor());
1064 
1065  // rows
1066  double check_delta_sum=0;
1067  for(i = 0; i <= ml_it_->size(); ++i) {
1068  pressure_diff(i, dofs_i, ele_type_i, delta_i, dirichlet_i);
1069  check_delta_sum+=delta_i;
1070  //columns
1071  for (j = 0; j <= ml_it_->size(); ++j) {
1072  pressure_diff(j, dofs_j, ele_type_j, delta_j, dirichlet_j);
1073 
1074  double scale = -master_sigma * delta_i * delta_j / delta_0;
1075  product = scale * tensor_average[ele_type_i][ele_type_j];
1076 
1077  arma::vec rhs(dofs_i.size());
1078  rhs.zeros();
1079  ls.set_values( dofs_i, dofs_j, product, rhs, dirichlet_i, dirichlet_j);
1080  //auto dofs_i_cp=dofs_i;
1081  //auto dofs_j_cp=dofs_j;
1082  //ls.set_values( dofs_i_cp, dofs_j_cp, product, rhs, dirichlet_i, dirichlet_j);
1083  }
1084  }
1085  OLD_ASSERT(check_delta_sum < 1E-5*delta_0, "sum err %f > 0\n", check_delta_sum/delta_0);
1086  } // loop over master elements
1087 }
1088 
1089 
1090 
1091  void P1_CouplingAssembler::add_sides(const Element * ele, unsigned int shift, vector<int> &dofs, vector<double> &dirichlet)
1092  {
1093 
1094  for(unsigned int i_side=0; i_side < ele->n_sides(); i_side++ ) {
1095  dofs[shift+i_side] = darcy_.mh_dh.row_4_edge[ele->side(i_side)->edge_idx()];
1096  Boundary * bcd = ele->side(i_side)->cond();
1097 
1098  if (bcd) {
1099  ElementAccessor<3> b_ele = bcd->element_accessor();
1100  auto type = (DarcyMH::EqData::BC_Type)darcy_.data_->bc_type.value(b_ele.centre(), b_ele);
1101  //DebugOut().fmt("bcd id: {} sidx: {} type: {}\n", ele->id(), i_side, type);
1102  if (type == DarcyMH::EqData::dirichlet) {
1103  //DebugOut().fmt("Dirichlet: {}\n", ele->index());
1104  dofs[shift + i_side] = -dofs[shift + i_side];
1105  double bc_pressure = darcy_.data_->bc_pressure.value(b_ele.centre(), b_ele);
1106  dirichlet[shift + i_side] = bc_pressure;
1107  }
1108  }
1109  }
1110  }
1111 
1112 
1113 /**
1114  * P1 connection of different dimensions
1115  *
1116  * - 20.11. 2014 - very poor convergence, big error in pressure even at internal part of the fracture
1117  */
1118 
1119 void P1_CouplingAssembler::assembly(LinSys &ls) {
1120 
1121  for (const Intersection &intersec : intersections_) {
1122  const Element * master = intersec.master_iter();
1123  const Element * slave = intersec.slave_iter();
1124 
1125  add_sides(slave, 0, dofs, dirichlet);
1126  add_sides(master, 3, dofs, dirichlet);
1127 
1128  double master_sigma=darcy_.data_->sigma.value( master->centre(), master->element_accessor());
1129 
1130 /*
1131  * Local coordinates on 1D
1132  * t0
1133  * node 0: 0.0
1134  * node 1: 1.0
1135  *
1136  * base fce points
1137  * t0 = 0.0 on side 0 node 0
1138  * t0 = 1.0 on side 1 node 1
1139  *
1140  * Local coordinates on 2D
1141  * t0 t1
1142  * node 0: 0.0 0.0
1143  * node 1: 1.0 0.0
1144  * node 2: 0.0 1.0
1145  *
1146  * base fce points
1147  * t0=0.5, t1=0.0 on side 0 nodes (0,1)
1148  * t0=0.5, t1=0.5 on side 1 nodes (1,2)
1149  * t0=0.0, t1=0.5 on side 2 nodes (2,0)
1150  */
1151 
1152 
1153 
1154  arma::vec point_Y(1);
1155  point_Y.fill(1.0);
1156  arma::vec point_2D_Y(intersec.map_to_slave(point_Y)); // local coordinates of Y on slave (1, t0, t1)
1157  arma::vec point_1D_Y(intersec.map_to_master(point_Y)); // local coordinates of Y on master (1, t0)
1158 
1159  arma::vec point_X(1);
1160  point_X.fill(0.0);
1161  arma::vec point_2D_X(intersec.map_to_slave(point_X)); // local coordinates of X on slave (1, t0, t1)
1162  arma::vec point_1D_X(intersec.map_to_master(point_X)); // local coordinates of X on master (1, t0)
1163 
1164  arma::mat base_2D(3, 3);
1165  // basis functions are numbered as sides
1166  // TODO:
1167  // Use RT finite element to evaluate following matrices.
1168 
1169  // Ravirat - Thomas base functions evaluated in points (0,0), (1,0), (0,1)
1170  // 2D RT_i(t0, t1) = a0 + a1*t0 + a2*t1
1171  // a0 a1 a2
1172  base_2D << 1.0 << 0.0 << -2.0 << arma::endr // RT for side 0
1173  << 1.0 << -2.0 << 0.0 << arma::endr // RT for side 1
1174  << -1.0 << 2.0 << 2.0 << arma::endr;// RT for side 2
1175 
1176 
1177  arma::mat base_1D(2, 2);
1178  // Ravirat - Thomas base functions evaluated in points (0,0), (1,0), (0,1)
1179  // 1D RT_i(t0) = a0 + a1 * t0
1180  // a0 a1
1181  base_1D << 1.0 << -1.0 << arma::endr // RT for side 0,
1182  << 0.0 << 1.0 << arma::endr; // RT for side 1,
1183 
1184 
1185 
1186  // Consider both 2D and 1D value are defined for the test function
1187  // related to the each of 5 DOFs involved in the intersection.
1188  // One of these values is always zero.
1189  // Compute difference of the 2D and 1D value for every DOF.
1190  // Compute value of this difference in both endpoints X,Y of the intersection.
1191 
1192  arma::vec difference_in_Y(5);
1193  arma::vec difference_in_X(5);
1194  // slave sides 0,1,2
1195  difference_in_Y.subvec(0, 2) = -base_2D * point_2D_Y;
1196  difference_in_X.subvec(0, 2) = -base_2D * point_2D_X;
1197  // master sides 3,4
1198  difference_in_Y.subvec(3, 4) = base_1D * point_1D_Y;
1199  difference_in_X.subvec(3, 4) = base_1D * point_1D_X;
1200 
1201  // applying the Simpson's rule
1202  // to the product of two linear functions f, g we get
1203  // (b-a)/6 * ( 3*f(a)*g(a) + 3*f(b)*g(b) + 2*f(a)*g(b) + 2*f(b)*g(a) )
1204  arma::mat A(5, 5);
1205  for (int i = 0; i < 5; ++i) {
1206  for (int j = 0; j < 5; ++j) {
1207  A(i, j) = -master_sigma * intersec.intersection_true_size() *
1208  ( difference_in_Y[i] * difference_in_Y[j]
1209  + difference_in_Y[i] * difference_in_X[j]/2
1210  + difference_in_X[i] * difference_in_Y[j]/2
1211  + difference_in_X[i] * difference_in_X[j]
1212  ) * (1.0 / 3.0);
1213 
1214  }
1215  }
1216  auto dofs_cp=dofs;
1217  ls.set_values( dofs_cp, dofs_cp, A, rhs, dirichlet, dirichlet);
1218 
1219  }
1220 }
1221 
1222 
1223 
1224 /*******************************************************************************
1225  * COMPOSE WATER MH MATRIX WITHOUT SCHUR COMPLEMENT
1226  ******************************************************************************/
1227 
1228 void DarcyMH::create_linear_system(Input::AbstractRecord in_rec) {
1229 
1230  START_TIMER("preallocation");
1231 
1232  if (schur0 == NULL) { // create Linear System for MH matrix
1233 
1234  if (in_rec.type() == LinSys_BDDC::get_input_type()) {
1235 #ifdef FLOW123D_HAVE_BDDCML
1236  WarningOut() << "For BDDC no Schur complements are used.";
1237  mh_dh.prepare_parallel_bddc();
1238  n_schur_compls = 0;
1239  LinSys_BDDC *ls = new LinSys_BDDC(mh_dh.global_row_4_sub_row->size(), &(*mh_dh.rows_ds),
1240  3, // 3 == la::BddcmlWrapper::SPD_VIA_SYMMETRICGENERAL
1241  1, // 1 == number of subdomains per process
1242  true); // swap signs of matrix and rhs to make the matrix SPD
1243  ls->set_from_input(in_rec);
1244  ls->set_solution( NULL );
1245  // possible initialization particular to BDDC
1246  START_TIMER("BDDC set mesh data");
1247  set_mesh_data_for_bddc(ls);
1248  schur0=ls;
1249  END_TIMER("BDDC set mesh data");
1250 #else
1251  xprintf(Err, "Flow123d was not build with BDDCML support.\n");
1252 #endif // FLOW123D_HAVE_BDDCML
1253  }
1254  else if (in_rec.type() == LinSys_PETSC::get_input_type()) {
1255  // use PETSC for serial case even when user wants BDDC
1256  if (n_schur_compls > 2) {
1257  WarningOut() << "Invalid number of Schur Complements. Using 2.";
1258  n_schur_compls = 2;
1259  }
1260 
1261  LinSys_PETSC *schur1, *schur2;
1262 
1263  if (n_schur_compls == 0) {
1264  LinSys_PETSC *ls = new LinSys_PETSC( &(*mh_dh.rows_ds) );
1265 
1266  // temporary solution; we have to set precision also for sequantial case of BDDC
1267  // final solution should be probably call of direct solver for oneproc case
1268  if (in_rec.type() != LinSys_BDDC::get_input_type()) ls->set_from_input(in_rec);
1269  else {
1270  ls->LinSys::set_from_input(in_rec); // get only common options
1271  }
1272 
1273  ls->set_solution( NULL );
1274  schur0=ls;
1275  } else {
1276  IS is;
1277  ISCreateStride(PETSC_COMM_WORLD, mh_dh.side_ds->lsize(), mh_dh.rows_ds->begin(), 1, &is);
1278  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
1279 
1280  SchurComplement *ls = new SchurComplement(is, &(*mh_dh.rows_ds));
1281 
1282  // make schur1
1283  Distribution *ds = ls->make_complement_distribution();
1284  if (n_schur_compls==1) {
1285  schur1 = new LinSys_PETSC(ds);
1286  schur1->set_positive_definite();
1287  } else {
1288  IS is;
1289  ISCreateStride(PETSC_COMM_WORLD, mh_dh.el_ds->lsize(), ls->get_distribution()->begin(), 1, &is);
1290  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
1291  SchurComplement *ls1 = new SchurComplement(is, ds); // is is deallocated by SchurComplement
1292  ls1->set_negative_definite();
1293 
1294  // make schur2
1295  schur2 = new LinSys_PETSC( ls1->make_complement_distribution() );
1296  schur2->set_positive_definite();
1297  ls1->set_complement( schur2 );
1298  schur1 = ls1;
1299  }
1300  ls->set_complement( schur1 );
1301  ls->set_from_input(in_rec);
1302  ls->set_solution( NULL );
1303  schur0=ls;
1304  }
1305 
1306  START_TIMER("PETSC PREALLOCATION");
1307  schur0->set_symmetric();
1308  schur0->start_allocation();
1309  assembly_mh_matrix(MultidimAssembler()); // preallocation
1310  VecZeroEntries(schur0->get_solution());
1311  END_TIMER("PETSC PREALLOCATION");
1312  }
1313  else {
1314  xprintf(Err, "Unknown solver type. Internal error.\n");
1315  }
1316 
1317  END_TIMER("preallocation");
1318  make_serial_scatter();
1319 
1320  }
1321 
1322 }
1323 
1324 
1325 
1326 
1327 void DarcyMH::assembly_linear_system() {
1328  START_TIMER("DarcyFlowMH_Steady::assembly_linear_system");
1329 
1330  is_linear_=true;
1331  bool is_steady = zero_time_term();
1332  //DebugOut() << "Assembly linear system\n";
1333  if (data_changed_) {
1334  data_changed_ = false;
1335  //DebugOut() << "Data changed\n";
1336  // currently we have no optimization for cases when just time term data or RHS data are changed
1337  START_TIMER("full assembly");
1338  if (typeid(*schur0) != typeid(LinSys_BDDC)) {
1339  schur0->start_add_assembly(); // finish allocation and create matrix
1340  }
1341 
1342  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(data_);
1343 
1344  schur0->mat_zero_entries();
1345  schur0->rhs_zero_entries();
1346 
1347  assembly_source_term();
1348  assembly_mh_matrix( multidim_assembler ); // fill matrix
1349 
1350  schur0->finish_assembly();
1351  schur0->set_matrix_changed();
1352  //MatView( *const_cast<Mat*>(schur0->get_matrix()), PETSC_VIEWER_STDOUT_WORLD );
1353  //VecView( *const_cast<Vec*>(schur0->get_rhs()), PETSC_VIEWER_STDOUT_WORLD);
1354 
1355  if (! is_steady) {
1356  START_TIMER("fix time term");
1357  //DebugOut() << "setup time term\n";
1358  // assembly time term and rhs
1359  setup_time_term();
1360  modify_system();
1361  }
1362  else if (balance_ != nullptr)
1363  {
1364  balance_->start_mass_assembly(water_balance_idx_);
1365  balance_->finish_mass_assembly(water_balance_idx_);
1366  }
1367  END_TIMER("full assembly");
1368  } else {
1369  START_TIMER("modify system");
1370  if (! is_steady) {
1371  modify_system();
1372  } else {
1373  //xprintf(PrgErr, "Planned computation time for steady solver, but data are not changed.\n");
1374  }
1375  END_TIMER("modiffy system");
1376  }
1377 
1378 }
1379 
1380 
1381 
1382 void DarcyMH::set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls) {
1383  START_TIMER("DarcyFlowMH_Steady::set_mesh_data_for_bddc");
1384  // prepare mesh for BDDCML
1385  // initialize arrays
1386  // auxiliary map for creating coordinates of local dofs and global-to-local numbering
1387  std::map<int,arma::vec3> localDofMap;
1388  // connectivity for the subdomain, i.e. global dof numbers on element, stored element-by-element
1389  // Indices of Nodes on Elements
1390  std::vector<int> inet;
1391  // number of degrees of freedom on elements - determines elementwise chunks of INET array
1392  // Number of Nodes on Elements
1393  std::vector<int> nnet;
1394  // Indices of Subdomain Elements in Global Numbering - for local elements, their global indices
1395  std::vector<int> isegn;
1396 
1397  // This array is currently not used in BDDCML, it was used as an interface scaling alternative to scaling
1398  // by diagonal. It corresponds to the rho-scaling.
1399  std::vector<double> element_permeability;
1400 
1401  // maximal and minimal dimension of elements
1402  uint elDimMax = 1;
1403  uint elDimMin = 3;
1404  for ( unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++ ) {
1405  auto ele_ac = mh_dh.accessor(i_loc);
1406  // for each element, create local numbering of dofs as fluxes (sides), pressure (element centre), Lagrange multipliers (edges), compatible connections
1407 
1408  elDimMax = std::max( elDimMax, ele_ac.dim() );
1409  elDimMin = std::min( elDimMin, ele_ac.dim() );
1410 
1411  isegn.push_back( ele_ac.ele_global_idx() );
1412  int nne = 0;
1413 
1414  FOR_ELEMENT_SIDES(ele_ac.full_iter(), si) {
1415  // insert local side dof
1416  int side_row = ele_ac.side_row(si);
1417  arma::vec3 coord = ele_ac.side(si)->centre();
1418 
1419  localDofMap.insert( std::make_pair( side_row, coord ) );
1420  inet.push_back( side_row );
1421  nne++;
1422  }
1423 
1424  // insert local pressure dof
1425  int el_row = ele_ac.ele_row();
1426  arma::vec3 coord = ele_ac.centre();
1427  localDofMap.insert( std::make_pair( el_row, coord ) );
1428  inet.push_back( el_row );
1429  nne++;
1430 
1431  FOR_ELEMENT_SIDES(ele_ac.full_iter(), si) {
1432  // insert local edge dof
1433  int edge_row = ele_ac.edge_row(si);
1434  arma::vec3 coord = ele_ac.side(si)->centre();
1435 
1436  localDofMap.insert( std::make_pair( edge_row, coord ) );
1437  inet.push_back( edge_row );
1438  nne++;
1439  }
1440 
1441  // insert dofs related to compatible connections
1442  for ( unsigned int i_neigh = 0; i_neigh < ele_ac.full_iter()->n_neighs_vb; i_neigh++) {
1443  int edge_row = mh_dh.row_4_edge[ ele_ac.full_iter()->neigh_vb[i_neigh]->edge_idx() ];
1444  arma::vec3 coord = ele_ac.full_iter()->neigh_vb[i_neigh]->edge()->side(0)->centre();
1445 
1446  localDofMap.insert( std::make_pair( edge_row, coord ) );
1447  inet.push_back( edge_row );
1448  nne++;
1449  }
1450 
1451  nnet.push_back( nne );
1452 
1453  // version for rho scaling
1454  // trace computation
1455  arma::vec3 centre = ele_ac.centre();
1456  double conduct = data_->conductivity.value( centre , ele_ac.element_accessor() );
1457  auto aniso = data_->anisotropy.value( centre, ele_ac.element_accessor() );
1458 
1459  // compute mean on the diagonal
1460  double coef = 0.;
1461  for ( int i = 0; i < 3; i++) {
1462  coef = coef + aniso.at(i,i);
1463  }
1464  // Maybe divide by cs
1465  coef = conduct*coef / 3;
1466 
1467  OLD_ASSERT( coef > 0.,
1468  "Zero coefficient of hydrodynamic resistance %f . \n ", coef );
1469  element_permeability.push_back( 1. / coef );
1470  }
1471  //convert set of dofs to vectors
1472  // number of nodes (= dofs) on the subdomain
1473  int numNodeSub = localDofMap.size();
1474  OLD_ASSERT_EQUAL( (unsigned int)numNodeSub, mh_dh.global_row_4_sub_row->size() );
1475  // Indices of Subdomain Nodes in Global Numbering - for local nodes, their global indices
1476  std::vector<int> isngn( numNodeSub );
1477  // pseudo-coordinates of local nodes (i.e. dofs)
1478  // they need not be exact, they are used just for some geometrical considerations in BDDCML,
1479  // such as selection of corners maximizing area of a triangle, bounding boxes fro subdomains to
1480  // find candidate neighbours etc.
1481  std::vector<double> xyz( numNodeSub * 3 ) ;
1482  int ind = 0;
1483  std::map<int,arma::vec3>::iterator itB = localDofMap.begin();
1484  for ( ; itB != localDofMap.end(); ++itB ) {
1485  isngn[ind] = itB -> first;
1486 
1487  arma::vec3 coord = itB -> second;
1488  for ( int j = 0; j < 3; j++ ) {
1489  xyz[ j*numNodeSub + ind ] = coord[j];
1490  }
1491 
1492  ind++;
1493  }
1494  localDofMap.clear();
1495 
1496  // Number of Nodal Degrees of Freedom
1497  // nndf is trivially one - dofs coincide with nodes
1498  std::vector<int> nndf( numNodeSub, 1 );
1499 
1500  // prepare auxiliary map for renumbering nodes
1501  typedef std::map<int,int> Global2LocalMap_; //! type for storage of global to local map
1502  Global2LocalMap_ global2LocalNodeMap;
1503  for ( unsigned ind = 0; ind < isngn.size(); ++ind ) {
1504  global2LocalNodeMap.insert( std::make_pair( static_cast<unsigned>( isngn[ind] ), ind ) );
1505  }
1506 
1507  // renumber nodes in the inet array to locals
1508  int indInet = 0;
1509  for ( unsigned int iEle = 0; iEle < isegn.size(); iEle++ ) {
1510  int nne = nnet[ iEle ];
1511  for ( int ien = 0; ien < nne; ien++ ) {
1512 
1513  int indGlob = inet[indInet];
1514  // map it to local node
1515  Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1516  OLD_ASSERT( pos != global2LocalNodeMap.end(),
1517  "Cannot remap node index %d to local indices. \n ", indGlob );
1518  int indLoc = static_cast<int> ( pos -> second );
1519 
1520  // store the node
1521  inet[ indInet++ ] = indLoc;
1522  }
1523  }
1524 
1525  int numNodes = size;
1526  int numDofsInt = size;
1527  int spaceDim = 3; // TODO: what is the proper value here?
1528  int meshDim = elDimMax;
1529 
1530  bddc_ls -> load_mesh( spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1531 }
1532 
1533 
1534 
1535 
1536 //=============================================================================
1537 // DESTROY WATER MH SYSTEM STRUCTURE
1538 //=============================================================================
1539 DarcyMH::~DarcyMH() {
1540  if (schur0 != NULL) {
1541  delete schur0;
1542  VecScatterDestroy(&par_to_all);
1543  }
1544 
1545  if (solution != NULL) {
1546  chkerr(VecDestroy(&sol_vec));
1547  delete [] solution;
1548  }
1549 
1550  if (output_object) delete output_object;
1551 
1552 
1553 
1554 }
1555 
1556 
1557 // ================================================
1558 // PARALLLEL PART
1559 //
1560 
1561 
1562 void DarcyMH::make_serial_scatter() {
1563  START_TIMER("prepare scatter");
1564  // prepare Scatter form parallel to sequantial in original numbering
1565  {
1566  IS is_loc;
1567  int i, *loc_idx; //, si;
1568 
1569  // create local solution vector
1570  solution = new double[size];
1571  VecCreateSeqWithArray(PETSC_COMM_SELF,1, size, solution,
1572  &(sol_vec));
1573 
1574  // create seq. IS to scatter par solutin to seq. vec. in original order
1575  // use essentialy row_4_id arrays
1576  loc_idx = new int [size];
1577  i = 0;
1578  FOR_ELEMENTS(mesh_, ele) {
1579  FOR_ELEMENT_SIDES(ele,si) {
1580  loc_idx[i++] = mh_dh.side_row_4_id[ mh_dh.side_dof( ele->side(si) ) ];
1581  }
1582  }
1583  FOR_ELEMENTS(mesh_, ele) {
1584  loc_idx[i++] = mh_dh.row_4_el[ele.index()];
1585  }
1586  for(unsigned int i_edg=0; i_edg < mesh_->n_edges(); i_edg++) {
1587  loc_idx[i++] = mh_dh.row_4_edge[i_edg];
1588  }
1589  OLD_ASSERT( i==size,"Size of array does not match number of fills.\n");
1590  //DBGPRINT_INT("loc_idx",size,loc_idx);
1591  ISCreateGeneral(PETSC_COMM_SELF, size, loc_idx, PETSC_COPY_VALUES, &(is_loc));
1592  delete [] loc_idx;
1593  VecScatterCreate(schur0->get_solution(), is_loc, sol_vec,
1594  PETSC_NULL, &par_to_all);
1595  ISDestroy(&(is_loc));
1596  }
1597  solution_changed_for_scatter=true;
1598 
1599  END_TIMER("prepare scatter");
1600 
1601 }
1602 
1603 
1604 /*
1605 void mat_count_off_proc_values(Mat m, Vec v) {
1606  int n, first, last;
1607  const PetscInt *cols;
1608  Distribution distr(v);
1609 
1610  int n_off = 0;
1611  int n_on = 0;
1612  int n_off_rows = 0;
1613  MatGetOwnershipRange(m, &first, &last);
1614  for (int row = first; row < last; row++) {
1615  MatGetRow(m, row, &n, &cols, PETSC_NULL);
1616  bool exists_off = false;
1617  for (int i = 0; i < n; i++)
1618  if (distr.get_proc(cols[i]) != distr.myp())
1619  n_off++, exists_off = true;
1620  else
1621  n_on++;
1622  if (exists_off)
1623  n_off_rows++;
1624  MatRestoreRow(m, row, &n, &cols, PETSC_NULL);
1625  }
1626 }
1627 */
1628 
1629 
1630 
1631 
1632 
1633 
1634 
1635 
1636 
1637 
1638 
1639 
1640 void DarcyMH::read_initial_condition()
1641 {
1642  double *local_sol = schur0->get_solution_array();
1643 
1644  // cycle over local element rows
1645 
1646  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1647  for (unsigned int i_loc_el = 0; i_loc_el < mh_dh.el_ds->lsize(); i_loc_el++) {
1648  auto ele_ac = mh_dh.accessor(i_loc_el);
1649  // set initial condition
1650  local_sol[ele_ac.ele_local_row()] = data_->init_pressure.value(ele_ac.centre(),ele_ac.element_accessor());
1651  }
1652 
1653  solution_changed_for_scatter=true;
1654 
1655 }
1656 
1657 void DarcyMH::setup_time_term() {
1658  // save diagonal of steady matrix
1659  MatGetDiagonal(*( schur0->get_matrix() ), steady_diagonal);
1660  // save RHS
1661  VecCopy(*( schur0->get_rhs()), steady_rhs);
1662 
1663 
1664  PetscScalar *local_diagonal;
1665  VecGetArray(new_diagonal,& local_diagonal);
1666 
1667  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1668 
1669  if (balance_ != nullptr)
1670  balance_->start_mass_assembly(water_balance_idx_);
1671 
1672  //DebugOut().fmt("time_term lsize: {} {}\n", mh_dh.el_ds->myp(), mh_dh.el_ds->lsize());
1673  for (unsigned int i_loc_el = 0; i_loc_el < mh_dh.el_ds->lsize(); i_loc_el++) {
1674  auto ele_ac = mh_dh.accessor(i_loc_el);
1675 
1676  // set new diagonal
1677  double diagonal_coeff = data_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor())
1678  * data_->storativity.value(ele_ac.centre(), ele_ac.element_accessor())
1679  * ele_ac.measure();
1680  local_diagonal[ele_ac.ele_local_row()]= - diagonal_coeff / time_->dt();
1681 
1682  //DebugOut().fmt("time_term: {} {} {} {} {}\n", mh_dh.el_ds->myp(), ele_ac.ele_global_idx(), i_loc_row, i_loc_el + mh_dh.side_ds->lsize(), diagonal_coeff);
1683  if (balance_ != nullptr)
1684  balance_->add_mass_matrix_values(water_balance_idx_,
1685  ele_ac.region().bulk_idx(), { int(ele_ac.ele_row()) }, {diagonal_coeff});
1686  }
1687  VecRestoreArray(new_diagonal,& local_diagonal);
1688  MatDiagonalSet(*( schur0->get_matrix() ), new_diagonal, ADD_VALUES);
1689 
1690  solution_changed_for_scatter=true;
1691  schur0->set_matrix_changed();
1692 
1693  if (balance_ != nullptr)
1694  balance_->finish_mass_assembly(water_balance_idx_);
1695 }
1696 
1697 void DarcyMH::modify_system() {
1698  START_TIMER("modify system");
1699  if (time_->is_changed_dt() && time_->step().index()>0) {
1700  double scale_factor=time_->step(-2).length()/time_->step().length();
1701  if (scale_factor != 1.0) {
1702  // if time step has changed and setup_time_term not called
1703  MatDiagonalSet(*( schur0->get_matrix() ),steady_diagonal, INSERT_VALUES);
1704 
1705  VecScale(new_diagonal, time_->last_dt()/time_->dt());
1706  MatDiagonalSet(*( schur0->get_matrix() ),new_diagonal, ADD_VALUES);
1707  schur0->set_matrix_changed();
1708  }
1709  }
1710 
1711  // modify RHS - add previous solution
1712  //VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, previous_solution);
1713  VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, schur0->get_solution());
1714  VecAXPY(*( schur0->get_rhs()), 1.0, steady_rhs);
1715  schur0->set_rhs_changed();
1716 
1717  //VecSwap(previous_solution, schur0->get_solution());
1718 }
1719 
1720 
1721 //-----------------------------------------------------------------------------
1722 // vim: set cindent:
TimeGovernor & time()
Definition: equation.hh:148
FieldSet & data()
Definition: equation.hh:193
void get_solution_vector(double *&vec, unsigned int &vec_size) override
FieldSet * eq_data_
Definition: equation.hh:232
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:42
void reinit(Mesh *mesh)
double measure() const
Definition: elements.cc:89
Output class for darcy_flow_mh model.
void assembly_mh_matrix(MultidimAssembler ma)
virtual void initialize_specific()
void get_parallel_solution_vector(Vec &vector) override
Accessor to input data conforming to declared Array.
Definition: accessors.hh:561
unsigned int uint
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:593
int is_linear_
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void make_intersec_elements()
Definition: mesh.cc:664
int tlevel() const
void init_eq_data()
static const Input::Type::Record & get_input_type()
The specification of output stream.
Definition: output_time.cc:38
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:56
Class for declaration of the input of type Bool.
Definition: type_base.hh:452
unsigned int edge_idx() const
Definition: side_impl.hh:61
virtual void assembly_source_term()
Source term is implemented differently in LMH version.
friend class DarcyFlowMHOutput
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:231
Boundary * cond() const
Definition: side_impl.hh:70
virtual void output_data() override
Write computed fields.
Wrappers for linear systems based on MPIAIJ and MATIS format.
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
bool is_end() const
Returns true if the actual time is greater than or equal to the end time.
void assembly(LinSys &ls)
void next_time()
Proceed to the next time according to current estimated time step.
static const int registrar
Registrar of class to factory.
unsigned int max_n_it_
VecScatter par_to_all
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:105
Vec steady_rhs
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
void initialize() override
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
Definition: mesh.h:95
std::vector< char > bc_switch_dirichlet
Idicator of dirichlet or neumann type of switch boundary conditions.
static const Input::Type::Record & type_field_descriptor()
Distribution * el_ds
virtual double get_solution_precision()=0
const RegionDB & region_db() const
Definition: mesh.h:147
void pressure_diff(int i_ele, vector< int > &dofs, unsigned int &ele_type, double &delta, arma::vec &dirichlet)
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
const TimeStep & step(int index=-1) const
static const std::string field_descriptor_record_description(const string &record_name)
Definition: field_common.cc:67
Class for declaration of the integral input data.
Definition: type_base.hh:483
double t() const
#define ADD_FIELD(name,...)
Definition: field_set.hh:269
Basic time management functionality for unsteady (and steady) solvers (class Equation).
unsigned int n_sides()
Definition: mesh.cc:170
bool solution_changed_for_scatter
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
std::shared_ptr< EqData > data_
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:301
Basic time management class.
virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it)=0
Class for declaration of inputs sequences.
Definition: type_base.hh:339
void view(const char *name="") const
MortarMethod mortar_method_
unsigned int min_n_it_
friend class P0_CouplingAssembler
unsigned int nonlinear_iteration_
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
void mat_set_value(int row, int col, double val)
Definition: linsys.hh:325
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:193
void update_solution() override
unsigned int size() const
Returns size of the container. This is independent of the allocated space.
Definition: sys_vector.hh:391
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:119
#define OLD_ASSERT(...)
Definition: global_defs.h:131
SideIter side()
Definition: boundaries.h:69
unsigned int n_elements() const
Definition: mesh.h:133
virtual void postprocess()
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
unsigned int dim() const
Definition: side_impl.hh:37
double * get_solution_array()
Definition: linsys.hh:281
Mesh & mesh()
Definition: equation.hh:176
std::shared_ptr< Balance > balance_
#define MPI_MIN
Definition: mpi.h:198
unsigned int edge_idx()
FLOW123D_FORCE_LINK_IN_CHILD(darcy_flow_mh)
int n_schur_compls
const Vec & get_solution()
Definition: linsys.hh:257
SideIter side()
Accessor to the data with type Type::Record.
Definition: accessors.hh:286
const Ret val(const string &key) const
unsigned int n_sides() const
#define xprintf(...)
Definition: system.hh:87
const Element * slave_iter() const
#define START_TIMER(tag)
Starts a timer with specified tag.
Mesh * mesh_
Definition: equation.hh:223
Element * element()
Definition: boundaries.cc:39
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
static const Input::Type::Instance & get_input_type()
bool data_changed_
double * solution
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:488
static Input::Type::Abstract & get_input_type()
ElementVector bc_elements
Definition: mesh.h:224
SideIter side(const unsigned int loc_index)
bool use_steady_assembly_
static const Input::Type::Record & get_input_type_specific()
static const Input::Type::Selection & get_bc_type_selection()
Return a Selection corresponding to enum BC_Type.
virtual void assembly_linear_system()
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:453
LinSys * schur0
unsigned int bc_ele_idx_
Definition: boundaries.h:76
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:213
DarcyMH(Mesh &mesh, const Input::Record in_rec)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
Definition: system.hh:59
Support classes for parallel programing.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
int set_upper_constraint(double upper, std::string message)
Sets upper constraint for the next time step estimating.
void set_values(int nrow, int *rows, int ncol, int *cols, PetscScalar *mat_vals, PetscScalar *rhs_vals)
Set values in the system matrix and values in the right-hand side vector on corresponding rows...
Definition: linsys.hh:342
void create_linear_system(Input::AbstractRecord rec)
Initialize global_row_4_sub_row.
MH_DofHandler mh_dh
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
const Selection & close() const
Close the Selection, no more values can be added.
double intersection_true_size() const
Input::Record input_record_
Definition: equation.hh:225
#define MPI_INT
Definition: mpi.h:160
double dt() const
std::shared_ptr< Distribution > rows_ds
LocalElementAccessorBase< 3 > accessor(uint local_ele_idx)
Vec steady_diagonal
#define ASSERT_DBG(expr)
Definition: asserts.hh:350
unsigned int water_balance_idx_
index of water balance within the Balance object.
Abstract linear system class.
static const Input::Type::Record & get_input_type()
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:234
BasicData Data
Definition: format.h:848
arma::vec::fixed< spacedim > centre() const
Definition: accessors.hh:91
#define OLD_ASSERT_EQUAL(a, b)
Definition: global_defs.h:133
EqData()
Creation of all fields.
mixed-hybrid model of linear Darcy flow, possibly unsteady.
void zero_time_step() override
virtual void read_initial_condition()
Record type proxy class.
Definition: type_record.hh:177
static Input::Type::Abstract & get_input_type()
Definition: linsys.cc:29
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
DarcyFlowMHOutput * output_object
unsigned int n_edges() const
Definition: mesh.h:141
virtual void prepare_new_time_step()
postprocess velocity field (add sources)
Class for representation SI units of Fields.
Definition: unit_si.hh:40
Vec new_diagonal
virtual double solution_precision() const
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
virtual double compute_residual()=0
Vec previous_solution
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:45
Template for classes storing finite set of named values.
friend class P1_CouplingAssembler
virtual int solve()=0
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:47
void rhs_set_value(int row, double val)
Definition: linsys.hh:337
double tolerance_
static const Input::Type::Record & get_input_type()
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
TimeGovernor * time_
Definition: equation.hh:224
virtual bool zero_time_term(bool time_global=false)
void output()
Calculate values for output.
unsigned int lsize(int proc) const
get local size