Flow123d  release_3.0.0-1159-g9e13a92
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/side_impl.hh"
37 #include "mesh/long_idx.hh"
38 #include "mesh/mesh.h"
39 #include "mesh/partitioning.hh"
40 #include "mesh/accessors.hh"
41 #include "mesh/range_wrapper.hh"
42 #include "la/distribution.hh"
43 #include "la/linsys.hh"
44 #include "la/linsys_PETSC.hh"
45 #include "la/linsys_BDDC.hh"
46 #include "la/schur.hh"
47 //#include "la/sparse_graph.hh"
49 
50 #include "flow/darcy_flow_mh.hh"
51 
53 
54 
55 #include "tools/time_governor.hh"
57 #include "fields/field.hh"
58 #include "fields/field_values.hh"
60 
61 #include "coupling/balance.hh"
62 
63 #include "la/vector_mpi.hh"
64 
65 #include "darcy_flow_assembly.hh"
66 
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", "No Mortar method is applied.")
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\n(zero normal flux over the boundary).")
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 (($ \\delta_d(q_d^N + \\sigma_d (h_d^R - h_d) )$)). "
97  "Specify the water inflow by the ``bc_flux`` field, the transition coefficient by ``bc_robin_sigma`` "
98  "and the reference pressure head or piezometric 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_d \\le h_d^D$)) and (($ -\\boldsymbol q_d\\cdot\\boldsymbol n \\le \\delta 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 the fields ``bc_switch_pressure`` (or ``bc_switch_piezo_head``) and ``bc_flux`` respectively."
106  )
107  .add_value(river, "river",
108  "River boundary condition. For the water level above the bedrock, (($H_d > H_d^S$)), the Robin boundary condition is used with the inflow given by: "
109  "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d) )$)). For the water level under the bedrock, constant infiltration is used: "
110  "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d^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", "Non-linear solver settings.")
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 solutions) to use.\nThis is usefull if the convergence criteria "
140  "does not characterize your goal well enough so it converges prematurely, possibly even without a single linear solution."
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 solutions) 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. By setting this flag to 'true', the solver "
147  "ends with convergence success on stagnation, but it reports warning about it.")
148  .close();
149 
150  return it::Record("Flow_Darcy_MH", "Mixed-Hybrid solver for 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("{}"),
157  "Non-linear solver for MH problem.")
158  .declare_key("output_stream", OutputTime::get_input_type(), it::Default("{}"),
159  "Output stream settings.\n Specify file format, precision etc.")
160 
161  .declare_key("output", DarcyFlowMHOutput::get_input_type(), IT::Default("{ \"fields\": [ \"pressure_p0\", \"velocity_p0\" ] }"),
162  "Specification of output fields and output times.")
164  "Output settings specific to Darcy flow model.\n"
165  "Includes raw output and some experimental functionality.")
166  .declare_key("balance", Balance::get_input_type(), it::Default("{}"),
167  "Settings for computing mass balance.")
169  "Time governor settings for the unsteady Darcy flow model.")
170  .declare_key("n_schurs", it::Integer(0,2), it::Default("2"),
171  "Number of Schur complements to perform when solving MH system.")
172  .declare_key("mortar_method", get_mh_mortar_selection(), it::Default("\"None\""),
173  "Method for coupling Darcy flow between dimensions on incompatible meshes. [Experimental]" )
174  .close();
175 }
176 
177 
178 const int DarcyMH::registrar =
179  Input::register_class< DarcyMH, Mesh &, const Input::Record >("Flow_Darcy_MH") +
181 
182 
183 
185 {
186  mortar_method_=NoMortar;
187 
188  *this += anisotropy.name("anisotropy")
189  .description("Anisotropy of the conductivity tensor.")
190  .input_default("1.0")
191  .units( UnitSI::dimensionless() );
192 
193  *this += cross_section.name("cross_section")
194  .description("Complement dimension parameter (cross section for 1D, thickness for 2D).")
195  .input_default("1.0")
196  .units( UnitSI().m(3).md() );
197 
198  *this += conductivity.name("conductivity")
199  .description("Isotropic conductivity scalar.")
200  .input_default("1.0")
201  .units( UnitSI().m().s(-1) )
202  .set_limits(0.0);
203 
204  *this += sigma.name("sigma")
205  .description("Transition coefficient between dimensions.")
206  .input_default("1.0")
207  .units( UnitSI::dimensionless() );
208 
209  *this += water_source_density.name("water_source_density")
210  .description("Water source density.")
211  .input_default("0.0")
212  .units( UnitSI().s(-1) );
213 
214  *this += bc_type.name("bc_type")
215  .description("Boundary condition type.")
216  .input_selection( get_bc_type_selection() )
217  .input_default("\"none\"")
218  .units( UnitSI::dimensionless() );
219 
220  *this += bc_pressure
221  .disable_where(bc_type, {none, seepage} )
222  .name("bc_pressure")
223  .description("Prescribed pressure value on the boundary. Used for all values of ``bc_type`` except ``none`` and ``seepage``. "
224  "See documentation of ``bc_type`` for exact meaning of ``bc_pressure`` in individual boundary condition types.")
225  .input_default("0.0")
226  .units( UnitSI().m() );
227 
228  *this += bc_flux
229  .disable_where(bc_type, {none, dirichlet} )
230  .name("bc_flux")
231  .description("Incoming water boundary flux. Used for bc_types : ``total_flux``, ``seepage``, ``river``.")
232  .input_default("0.0")
233  .units( UnitSI().m().s(-1) );
234 
235  *this += bc_robin_sigma
236  .disable_where(bc_type, {none, dirichlet, seepage} )
237  .name("bc_robin_sigma")
238  .description("Conductivity coefficient in the ``total_flux`` or the ``river`` boundary condition type.")
239  .input_default("0.0")
240  .units( UnitSI().s(-1) );
241 
242  *this += bc_switch_pressure
243  .disable_where(bc_type, {none, dirichlet, total_flux} )
244  .name("bc_switch_pressure")
245  .description("Critical switch pressure for ``seepage`` and ``river`` boundary conditions.")
246  .input_default("0.0")
247  .units( UnitSI().m() );
248 
249 
250  //these are for unsteady
251  *this += init_pressure.name("init_pressure")
252  .description("Initial condition for pressure in time dependent problems.")
253  .input_default("0.0")
254  .units( UnitSI().m() );
255 
256  *this += storativity.name("storativity")
257  .description("Storativity (in time dependent problems).")
258  .input_default("0.0")
259  .units( UnitSI().m(-1) );
260 
261  //time_term_fields = this->subset({"storativity"});
262  //main_matrix_fields = this->subset({"anisotropy", "conductivity", "cross_section", "sigma", "bc_type", "bc_robin_sigma"});
263  //rhs_fields = this->subset({"water_source_density", "bc_pressure", "bc_flux"});
264 }
265 
266 
267 
268 
269 
270 
271 //=============================================================================
272 // CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
273 // - do it in parallel:
274 // - initial distribution of elements, edges
275 //
276 /*! @brief CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
277  *
278  * Parameters {Solver,NSchurs} number of performed Schur
279  * complements (0,1,2) for water flow MH-system
280  *
281  */
282 //=============================================================================
283 DarcyMH::DarcyMH(Mesh &mesh_in, const Input::Record in_rec)
284 : DarcyFlowInterface(mesh_in, in_rec),
285  output_object(nullptr),
286  solution(nullptr),
287  data_changed_(false),
288  schur0(nullptr),
289  steady_diagonal(nullptr),
290  steady_rhs(nullptr),
291  new_diagonal(nullptr),
292  previous_solution(nullptr)
293 {
294 
295  START_TIMER("Darcy constructor");
296  {
297  auto time_record = input_record_.val<Input::Record>("time");
298  //if ( in_rec.opt_val("time", time_record) )
299  time_ = new TimeGovernor(time_record);
300  //else
301  // time_ = new TimeGovernor();
302  }
303 
304  data_ = make_shared<EqData>();
306 
307  data_->is_linear=true;
308 
309  size = mesh_->n_elements() + mesh_->n_sides() + mesh_->n_edges();
310  n_schur_compls = in_rec.val<int>("n_schurs");
311  data_->mortar_method_= in_rec.val<MortarMethod>("mortar_method");
312  if (data_->mortar_method_ != NoMortar) {
314  }
315 
316 
317 
318  //side_ds->view( std::cout );
319  //el_ds->view( std::cout );
320  //edge_ds->view( std::cout );
321  //rows_ds->view( std::cout );
322 
323 }
324 
325 
326 
328 //connecting data fields with mesh
329 {
330 
331  START_TIMER("data init");
332  data_->mesh = mesh_;
333  data_->mh_dh = &mh_dh;
334  data_->set_mesh(*mesh_);
335 
336  auto gravity_array = input_record_.val<Input::Array>("gravity");
337  std::vector<double> gvec;
338  gravity_array.copy_to(gvec);
339  gvec.push_back(0.0); // zero pressure shift
340  data_->gravity_ = arma::vec(gvec);
341  data_->gravity_vec_ = data_->gravity_.subvec(0,2);
342 
343  data_->bc_pressure.add_factory(
344  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
345  (data_->gravity_, "bc_piezo_head") );
346  data_->bc_switch_pressure.add_factory(
347  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
348  (data_->gravity_, "bc_switch_piezo_head") );
349  data_->init_pressure.add_factory(
350  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
351  (data_->gravity_, "init_piezo_head") );
352 
353 
354  data_->set_input_list( this->input_record_.val<Input::Array>("input_fields"), *time_ );
355  // Check that the time step was set for the transient simulation.
356  if (! zero_time_term(true) && time_->is_default() ) {
357  //THROW(ExcAssertMsg());
358  //THROW(ExcMissingTimeGovernor() << input_record_.ei_address());
359  MessageOut() << "Missing the key 'time', obligatory for the transient problems." << endl;
360  ASSERT(false);
361  }
362 
363  data_->mark_input_times(*time_);
364 }
365 
366 
367 
369 
370  init_eq_data();
372 
373  mh_dh.reinit(mesh_);
374  // Initialize bc_switch_dirichlet to size of global boundary.
375  data_->bc_switch_dirichlet.resize(mesh_->n_elements()+mesh_->n_elements(true), 1);
376 
377 
380  .val<Input::Record>("nonlinear_solver")
381  .val<Input::AbstractRecord>("linear_solver");
382 
383  // auxiliary set_time call since allocation assembly evaluates fields as well
386 
387 
388 
389  // allocate time term vectors
390  VecDuplicate(schur0->get_solution(), &previous_solution);
391  VecCreateMPI(PETSC_COMM_WORLD, mh_dh.rows_ds->lsize(),PETSC_DETERMINE,&(steady_diagonal));
392  VecDuplicate(steady_diagonal,& new_diagonal);
393  VecZeroEntries(new_diagonal);
394  VecDuplicate(steady_diagonal, &steady_rhs);
395 
396 
397  // initialization of balance object
398  balance_ = std::make_shared<Balance>("water", mesh_);
399  balance_->init_from_input(input_record_.val<Input::Record>("balance"), time());
400  data_->water_balance_idx = balance_->add_quantity("water_volume");
401  balance_->allocate(mh_dh.rows_ds->lsize(), 1);
402  balance_->units(UnitSI().m(3));
403 
404 
405  data_->balance = balance_;
406  data_->lin_sys = schur0;
407 
408 
410 }
411 
413 {}
414 
416 {
417 
418  /* TODO:
419  * - Allow solution reconstruction (pressure and velocity) from initial condition on user request.
420  * - Steady solution as an intitial condition may be forced by setting inti_time =-1, and set data for the steady solver in that time.
421  * Solver should be able to switch from and to steady case depending on the zero time term.
422  */
423 
425 
426  // zero_time_term means steady case
427  bool zero_time_term_from_right = zero_time_term();
428 
429 
430  if (zero_time_term_from_right) {
431  // steady case
432  VecZeroEntries(schur0->get_solution());
433  //read_initial_condition(); // Possible solution guess for steady case.
434  use_steady_assembly_ = true;
435  solve_nonlinear(); // with right limit data
436  } else {
437  VecZeroEntries(schur0->get_solution());
438  VecZeroEntries(previous_solution);
440  assembly_linear_system(); // in particular due to balance
441 // print_matlab_matrix("matrix_zero");
442  // TODO: reconstruction of solution in zero time.
443  }
444  //solution_output(T,right_limit); // data for time T in any case
445  output_data();
446 }
447 
448 //=============================================================================
449 // COMPOSE and SOLVE WATER MH System possibly through Schur complements
450 //=============================================================================
452 {
453  START_TIMER("Solving MH system");
454 
455  time_->next_time();
456 
457  time_->view("DARCY"); //time governor information output
459  bool zero_time_term_from_left=zero_time_term();
460 
461  bool jump_time = data_->storativity.is_jump_time();
462  if (! zero_time_term_from_left) {
463  // time term not treated as zero
464  // Unsteady solution up to the T.
465 
466  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
467  use_steady_assembly_ = false;
468  prepare_new_time_step(); //SWAP
469 
470  solve_nonlinear(); // with left limit data
471  if (jump_time) {
472  WarningOut() << "Output of solution discontinuous in time not supported yet.\n";
473  //solution_output(T, left_limit); // output use time T- delta*dt
474  //output_data();
475  }
476  }
477 
478  if (time_->is_end()) {
479  // output for unsteady case, end_time should not be the jump time
480  // but rether check that
481  if (! zero_time_term_from_left && ! jump_time) output_data();
482  return;
483  }
484 
486  bool zero_time_term_from_right=zero_time_term();
487  if (zero_time_term_from_right) {
488  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
489  use_steady_assembly_ = true;
490  solve_nonlinear(); // with right limit data
491 
492  } else if (! zero_time_term_from_left && jump_time) {
493  WarningOut() << "Discontinuous time term not supported yet.\n";
494  //solution_transfer(); // internally call set_time(T, left) and set_time(T,right) again
495  //solve_nonlinear(); // with right limit data
496  }
497  //solution_output(T,right_limit); // data for time T in any case
498  output_data();
499 
500 }
501 
502 bool DarcyMH::zero_time_term(bool time_global) {
503  if (time_global) {
504  return (data_->storativity.input_list_size() == 0);
505  } else {
506  return data_->storativity.field_result(mesh_->region_db().get_region_set("BULK")) == result_zeros;
507  }
508 }
509 
510 
512 {
513 
515  double residual_norm = schur0->compute_residual();
517  MessageOut().fmt("[nonlinear solver] norm of initial residual: {}\n", residual_norm);
518 
519  // Reduce is_linear flag.
520  int is_linear_common;
521  MPI_Allreduce(&(data_->is_linear), &is_linear_common,1, MPI_INT ,MPI_MIN,PETSC_COMM_WORLD);
522 
523  Input::Record nl_solver_rec = input_record_.val<Input::Record>("nonlinear_solver");
524  this->tolerance_ = nl_solver_rec.val<double>("tolerance");
525  this->max_n_it_ = nl_solver_rec.val<unsigned int>("max_it");
526  this->min_n_it_ = nl_solver_rec.val<unsigned int>("min_it");
527  if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->max_n_it_;
528 
529  if (! is_linear_common) {
530  // set tolerances of the linear solver unless they are set by user.
531  schur0->set_tolerances(0.1*this->tolerance_, 0.01*this->tolerance_, 100);
532  }
533  vector<double> convergence_history;
534 
535  Vec save_solution;
536  VecDuplicate(schur0->get_solution(), &save_solution);
537  while (nonlinear_iteration_ < this->min_n_it_ ||
538  (residual_norm > this->tolerance_ && nonlinear_iteration_ < this->max_n_it_ )) {
539  OLD_ASSERT_EQUAL( convergence_history.size(), nonlinear_iteration_ );
540  convergence_history.push_back(residual_norm);
541 
542  // stagnation test
543  if (convergence_history.size() >= 5 &&
544  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
545  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
546  // stagnation
547  if (input_record_.val<Input::Record>("nonlinear_solver").val<bool>("converge_on_stagnation")) {
548  WarningOut().fmt("Accept solution on stagnation. Its: {} Residual: {}\n", nonlinear_iteration_, residual_norm);
549  break;
550  } else {
551  THROW(ExcSolverDiverge() << EI_Reason("Stagnation."));
552  }
553  }
554 
555  if (! is_linear_common)
556  VecCopy( schur0->get_solution(), save_solution);
559 
560  // hack to make BDDC work with empty compute_residual
561  if (is_linear_common){
562  // we want to print this info in linear (and steady) case
563  residual_norm = schur0->compute_residual();
564  MessageOut().fmt("[nonlinear solver] lin. it: {}, reason: {}, residual: {}\n",
565  si.n_iterations, si.converged_reason, residual_norm);
566  break;
567  }
568  data_changed_=true; // force reassembly for non-linear case
569 
570  double alpha = 1; // how much of new solution
571  VecAXPBY(schur0->get_solution(), (1-alpha), alpha, save_solution);
572 
573  /*
574  double * sol;
575  unsigned int sol_size;
576  get_solution_vector(sol, sol_size);
577  if (mh_dh.el_ds->myp() == 0)
578  VecView(sol_vec, PETSC_VIEWER_STDOUT_SELF);
579  */
580 
581  //LogOut().fmt("Linear solver ended with reason: {} \n", si.converged_reason );
582  //OLD_ASSERT( si.converged_reason >= 0, "Linear solver failed to converge. Convergence reason %d \n", si.converged_reason );
584 
585  residual_norm = schur0->compute_residual();
586  MessageOut().fmt("[nonlinear solver] it: {} lin. it: {}, reason: {}, residual: {}\n",
587  nonlinear_iteration_, si.n_iterations, si.converged_reason, residual_norm);
588  }
589  chkerr(VecDestroy(&save_solution));
590  this -> postprocess();
591 
592  // adapt timestep
593  if (! this->zero_time_term()) {
594  double mult = 1.0;
595  if (nonlinear_iteration_ < 3) mult = 1.6;
596  if (nonlinear_iteration_ > 7) mult = 0.7;
597  int result = time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
598  //DebugOut().fmt("time adaptivity, res: {} it: {} m: {} dt: {} edt: {}\n", result, nonlinear_iteration_, mult, time_->dt(), time_->estimate_dt());
599  }
600 
602 
603 }
604 
605 
607 {
608  //VecSwap(previous_solution, schur0->get_solution());
609 }
610 
612 {
613  START_TIMER("postprocess");
614 
615  //fix velocity when mortar method is used
616  if(data_->mortar_method_ != MortarMethod::NoMortar){
617  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(data_);
618  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
619  auto ele_ac = mh_dh.accessor(i_loc);
620  unsigned int dim = ele_ac.dim();
621  multidim_assembler[dim-1]->fix_velocity(ele_ac);
622  }
623  }
624  //ElementAccessor<3> ele;
625 
626  // modify side fluxes in parallel
627  // for every local edge take time term on digonal and add it to the corresponding flux
628  /*
629  for (unsigned int i_loc = 0; i_loc < el_ds->lsize(); i_loc++) {
630  ele = mesh_->element_accessor(el_4_loc[i_loc]);
631  for (unsigned int i=0; i<ele->n_sides(); i++) {
632  side_rows[i] = side_row_4_id[ mh_dh.side_dof( ele_ac.side(i) ) ];
633  values[i] = -1.0 * ele_ac.measure() *
634  data.cross_section.value(ele_ac.centre(), ele_ac.element_accessor()) *
635  data.water_source_density.value(ele_ac.centre(), ele_ac.element_accessor()) /
636  ele_ac.n_sides();
637  }
638  VecSetValues(schur0->get_solution(), ele_ac.n_sides(), side_rows, values, ADD_VALUES);
639  }
640  VecAssemblyBegin(schur0->get_solution());
641  VecAssemblyEnd(schur0->get_solution());
642  */
643 }
644 
645 
647  START_TIMER("Darcy output data");
648  //time_->view("DARCY"); //time governor information output
649  this->output_object->output();
650 
651 
652  START_TIMER("Darcy balance output");
653  balance_->calculate_cumulative(data_->water_balance_idx, schur0->get_solution());
654  balance_->calculate_instant(data_->water_balance_idx, schur0->get_solution());
655  balance_->output();
656 }
657 
658 
660 {
661  return schur0->get_solution_precision();
662 }
663 
664 
665 
666 void DarcyMH::get_solution_vector(double * &vec, unsigned int &vec_size)
667 {
668  // TODO: make class for vectors (wrapper for PETSC or other) derived from LazyDependency
669  // and use its mechanism to manage dependency between vectors
671 
672  // scatter solution to all procs
673  VecScatterBegin(par_to_all, schur0->get_solution(), sol_vec, INSERT_VALUES, SCATTER_FORWARD);
674  VecScatterEnd( par_to_all, schur0->get_solution(), sol_vec, INSERT_VALUES, SCATTER_FORWARD);
676  }
677 
678  vec = solution;
679  vec_size = this->size;
680  OLD_ASSERT(vec != NULL, "Requested solution is not allocated!\n");
681 }
682 
684 {
685  vec=schur0->get_solution();
686  OLD_ASSERT(vec != NULL, "Requested solution is not allocated!\n");
687 }
688 
689 
690 // ===========================================================================================
691 //
692 // MATRIX ASSEMBLY - we use abstract assembly routine, where LS Mat/Vec SetValues
693 // are in fact pointers to allocating or filling functions - this is governed by Linsystem roitunes
694 //
695 // =======================================================================================
697 {
698  START_TIMER("DarcyFlowMH_Steady::assembly_steady_mh_matrix");
699 
700  // set auxiliary flag for switchting Dirichlet like BC
701  data_->force_bc_switch = use_steady_assembly_ && (nonlinear_iteration_ == 0);
702  data_->n_schur_compls = n_schur_compls;
703 
704 
705  balance_->start_flux_assembly(data_->water_balance_idx);
706 
707  // TODO: try to move this into balance, or have it in the generic assembler class, that should perform the cell loop
708  // including various pre- and post-actions
709  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
710  auto ele_ac = mh_dh.accessor(i_loc);
711  unsigned int dim = ele_ac.dim();
712  assembler[dim-1]->assemble(ele_ac);
713  }
714 
715 
716  balance_->finish_flux_assembly(data_->water_balance_idx);
717 
718 }
719 
720 
722 {
723  START_TIMER("DarcyFlowMH_Steady::allocate_mh_matrix");
724 
725  // set auxiliary flag for switchting Dirichlet like BC
726  data_->n_schur_compls = n_schur_compls;
727  LinSys *ls = schur0;
728 
729 
730 
731  int local_dofs[10];
732 
733  // to make space for second schur complement, max. 10 neighbour edges of one el.
734  double zeros[100000];
735  for(int i=0; i<100000; i++) zeros[i] = 0.0;
736 
737  std::vector<int> tmp_rows;
738  tmp_rows.reserve(200);
739 
740  unsigned int nsides, loc_size;
741 
742  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
743  auto ele_ac = mh_dh.accessor(i_loc);
744  nsides = ele_ac.n_sides();
745 
746  //allocate at once matrix [sides,ele,edges]x[sides,ele,edges]
747  loc_size = 1 + 2*nsides;
748  unsigned int i_side = 0;
749 
750  for (; i_side < nsides; i_side++) {
751  local_dofs[i_side] = ele_ac.side_row(i_side);
752  local_dofs[i_side+nsides] = ele_ac.edge_row(i_side);
753  }
754  local_dofs[i_side+nsides] = ele_ac.ele_row();
755  int * edge_rows = local_dofs + nsides;
756  //int ele_row = local_dofs[0];
757 
758  // whole local MH matrix
759  ls->mat_set_values(loc_size, local_dofs, loc_size, local_dofs, zeros);
760 
761 
762  // compatible neighborings rows
763  unsigned int n_neighs = ele_ac.element_accessor()->n_neighs_vb();
764  for (unsigned int i = 0; i < n_neighs; i++) {
765  // every compatible connection adds a 2x2 matrix involving
766  // current element pressure and a connected edge pressure
767  Neighbour *ngh = ele_ac.element_accessor()->neigh_vb[i];
768  int neigh_edge_row = mh_dh.row_4_edge[ ngh->edge_idx() ];
769  tmp_rows.push_back(neigh_edge_row);
770  //DebugOut() << "CC" << print_var(tmp_rows[i]);
771  }
772 
773  // allocate always also for schur 2
774  ls->mat_set_values(nsides+1, edge_rows, n_neighs, tmp_rows.data(), zeros); // (edges, ele) x (neigh edges)
775  ls->mat_set_values(n_neighs, tmp_rows.data(), nsides+1, edge_rows, zeros); // (neigh edges) x (edges, ele)
776  ls->mat_set_values(n_neighs, tmp_rows.data(), n_neighs, tmp_rows.data(), zeros); // (neigh edges) x (neigh edges)
777 
778  tmp_rows.clear();
779 
780  if (data_->mortar_method_ != NoMortar) {
781  auto &isec_list = mesh_->mixed_intersections().element_intersections_[ele_ac.ele_global_idx()];
782  for(auto &isec : isec_list ) {
783  IntersectionLocalBase *local = isec.second;
784  ElementAccessor<3> slave_ele = mesh_->element_accessor( local->bulk_ele_idx() );
785  //DebugOut().fmt("Alloc: {} {}", ele_ac.ele_global_idx(), local->bulk_ele_idx());
786  for(unsigned int i_side=0; i_side < slave_ele->n_sides(); i_side++) {
787  tmp_rows.push_back( mh_dh.row_4_edge[ slave_ele.side(i_side)->edge_idx() ] );
788  //DebugOut() << "aedge" << print_var(tmp_rows[i_rows-1]);
789  }
790  }
791  }
792  /*
793  for(unsigned int i_side=0; i_side < ele_ac.element_accessor()->n_sides(); i_side++) {
794  DebugOut() << "aedge:" << print_var(edge_rows[i_side]);
795  }*/
796 
797  ls->mat_set_values(nsides, edge_rows, tmp_rows.size(), tmp_rows.data(), zeros); // master edges x neigh edges
798  ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), nsides, edge_rows, zeros); // neigh edges x master edges
799  ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), tmp_rows.size(), tmp_rows.data(), zeros); // neigh edges x neigh edges
800 
801  }
802 /*
803  // alloc edge diagonal entries
804  if(rank == 0)
805  for( vector<Edge>::iterator edg = mesh_->edges.begin(); edg != mesh_->edges.end(); ++edg) {
806  int edg_idx = mh_dh.row_4_edge[edg->side(0)->edge_idx()];
807 
808 // for( vector<Edge>::iterator edg2 = mesh_->edges.begin(); edg2 != mesh_->edges.end(); ++edg2){
809 // int edg_idx2 = mh_dh.row_4_edge[edg2->side(0)->edge_idx()];
810 // if(edg_idx == edg_idx2){
811 // DBGCOUT(<< "P[ " << rank << " ] " << "edg alloc: " << edg_idx << " " << edg_idx2 << "\n");
812  ls->mat_set_value(edg_idx, edg_idx, 0.0);
813 // }
814 // }
815  }
816  */
817  /*
818  if (mortar_method_ == MortarP0) {
819  P0_CouplingAssembler(*this).assembly(*ls);
820  } else if (mortar_method_ == MortarP1) {
821  P1_CouplingAssembler(*this).assembly(*ls);
822  }*/
823 }
824 
826 {
827  START_TIMER("assembly source term");
828  balance_->start_source_assembly(data_->water_balance_idx);
829 
830  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
831  auto ele_ac = mh_dh.accessor(i_loc);
832 
833  double cs = data_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor());
834 
835  // set sources
836  double source = ele_ac.measure() * cs *
837  data_->water_source_density.value(ele_ac.centre(), ele_ac.element_accessor());
838  schur0->rhs_set_value(ele_ac.ele_row(), -1.0 * source );
839 
840  balance_->add_source_values(data_->water_balance_idx, ele_ac.region().bulk_idx(), {(LongIdx) ele_ac.ele_local_row()}, {0}, {source});
841  }
842 
843  balance_->finish_source_assembly(data_->water_balance_idx);
844 }
845 
846 
847 
848 
849 /*******************************************************************************
850  * COMPOSE WATER MH MATRIX WITHOUT SCHUR COMPLEMENT
851  ******************************************************************************/
852 
854 
855  START_TIMER("preallocation");
856 
857  if (schur0 == NULL) { // create Linear System for MH matrix
858 
859  if (in_rec.type() == LinSys_BDDC::get_input_type()) {
860 #ifdef FLOW123D_HAVE_BDDCML
861  WarningOut() << "For BDDC no Schur complements are used.";
863  n_schur_compls = 0;
865  3, // 3 == la::BddcmlWrapper::SPD_VIA_SYMMETRICGENERAL
866  1, // 1 == number of subdomains per process
867  true); // swap signs of matrix and rhs to make the matrix SPD
868  ls->set_from_input(in_rec);
869  ls->set_solution();
870  // possible initialization particular to BDDC
871  START_TIMER("BDDC set mesh data");
873  schur0=ls;
874  END_TIMER("BDDC set mesh data");
875 #else
876  xprintf(Err, "Flow123d was not build with BDDCML support.\n");
877 #endif // FLOW123D_HAVE_BDDCML
878  }
879  else if (in_rec.type() == LinSys_PETSC::get_input_type()) {
880  // use PETSC for serial case even when user wants BDDC
881  if (n_schur_compls > 2) {
882  WarningOut() << "Invalid number of Schur Complements. Using 2.";
883  n_schur_compls = 2;
884  }
885 
886  LinSys_PETSC *schur1, *schur2;
887 
888  if (n_schur_compls == 0) {
889  LinSys_PETSC *ls = new LinSys_PETSC( &(*mh_dh.rows_ds) );
890 
891  // temporary solution; we have to set precision also for sequantial case of BDDC
892  // final solution should be probably call of direct solver for oneproc case
893  if (in_rec.type() != LinSys_BDDC::get_input_type()) ls->set_from_input(in_rec);
894  else {
895  ls->LinSys::set_from_input(in_rec); // get only common options
896  }
897 
898  ls->set_solution();
899  schur0=ls;
900  } else {
901  IS is;
902  ISCreateStride(PETSC_COMM_WORLD, mh_dh.side_ds->lsize(), mh_dh.rows_ds->begin(), 1, &is);
903  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
904 
905  SchurComplement *ls = new SchurComplement(&(*mh_dh.rows_ds), is);
906 
907  // make schur1
908  Distribution *ds = ls->make_complement_distribution();
909  if (n_schur_compls==1) {
910  schur1 = new LinSys_PETSC(ds);
911  schur1->set_positive_definite();
912  } else {
913  IS is;
914  ISCreateStride(PETSC_COMM_WORLD, mh_dh.el_ds->lsize(), ls->get_distribution()->begin(), 1, &is);
915  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
916  SchurComplement *ls1 = new SchurComplement(ds, is); // is is deallocated by SchurComplement
917  ls1->set_negative_definite();
918 
919  // make schur2
920  schur2 = new LinSys_PETSC( ls1->make_complement_distribution() );
921  schur2->set_positive_definite();
922  ls1->set_complement( schur2 );
923  schur1 = ls1;
924  }
925  ls->set_complement( schur1 );
926  ls->set_from_input(in_rec);
927  ls->set_solution();
928  schur0=ls;
929  }
930 
931  START_TIMER("PETSC PREALLOCATION");
934 
936 
937  VecZeroEntries(schur0->get_solution());
938  END_TIMER("PETSC PREALLOCATION");
939  }
940  else {
941  xprintf(Err, "Unknown solver type. Internal error.\n");
942  }
943 
944  END_TIMER("preallocation");
946 
947  }
948 
949 }
950 
951 
952 
953 
955  START_TIMER("DarcyFlowMH_Steady::assembly_linear_system");
956 
957  data_->is_linear=true;
958  bool is_steady = zero_time_term();
959  //DebugOut() << "Assembly linear system\n";
960  if (data_changed_) {
961  data_changed_ = false;
962  //DebugOut() << "Data changed\n";
963  // currently we have no optimization for cases when just time term data or RHS data are changed
964  START_TIMER("full assembly");
965  if (typeid(*schur0) != typeid(LinSys_BDDC)) {
966  schur0->start_add_assembly(); // finish allocation and create matrix
967  }
968 
971 
973 
974  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(data_);
975  assembly_mh_matrix( multidim_assembler ); // fill matrix
976 
978 // print_matlab_matrix("matrix");
980  //MatView( *const_cast<Mat*>(schur0->get_matrix()), PETSC_VIEWER_STDOUT_WORLD );
981  //VecView( *const_cast<Vec*>(schur0->get_rhs()), PETSC_VIEWER_STDOUT_WORLD);
982 
983  if (! is_steady) {
984  START_TIMER("fix time term");
985  //DebugOut() << "setup time term\n";
986  // assembly time term and rhs
987  setup_time_term();
988  modify_system();
989  }
990  else
991  {
992  balance_->start_mass_assembly(data_->water_balance_idx);
993  balance_->finish_mass_assembly(data_->water_balance_idx);
994  }
995  END_TIMER("full assembly");
996  } else {
997  START_TIMER("modify system");
998  if (! is_steady) {
999  modify_system();
1000  } else {
1001  //xprintf(PrgErr, "Planned computation time for steady solver, but data are not changed.\n");
1002  }
1003  END_TIMER("modiffy system");
1004  }
1005 
1006 }
1007 
1008 
1009 void DarcyMH::print_matlab_matrix(std::string matlab_file)
1010 {
1011  std::string output_file;
1012 
1013  if ( typeid(*schur0) == typeid(LinSys_BDDC) ){
1014 // WarningOut() << "Can output matrix only on a single processor.";
1015 // output_file = FilePath(matlab_file + "_bddc.m", FilePath::output_file);
1016 // ofstream os( output_file );
1017 // auto bddc = static_cast<LinSys_BDDC*>(schur0);
1018 // bddc->print_matrix(os);
1019  }
1020  else {//if ( typeid(*schur0) == typeid(LinSys_PETSC) ){
1021  output_file = FilePath(matlab_file + ".m", FilePath::output_file);
1022  PetscViewer viewer;
1023  PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
1024  PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1025  MatView( *const_cast<Mat*>(schur0->get_matrix()), viewer);
1026  VecView( *const_cast<Vec*>(schur0->get_rhs()), viewer);
1027  }
1028 // else{
1029 // WarningOut() << "No matrix output available for the current solver.";
1030 // return;
1031 // }
1032 
1033  // compute h_min for different dimensions
1034  double d_max = std::numeric_limits<double>::max();
1035  double h1 = d_max, h2 = d_max, h3 = d_max;
1036  double he2 = d_max, he3 = d_max;
1037  for (auto ele : mesh_->elements_range()) {
1038  switch(ele->dim()){
1039  case 1: h1 = std::min(h1,ele.measure()); break;
1040  case 2: h2 = std::min(h2,ele.measure()); break;
1041  case 3: h3 = std::min(h3,ele.measure()); break;
1042  }
1043 
1044  for (unsigned int j=0; j<ele->n_sides(); j++) {
1045  switch(ele->dim()){
1046  case 2: he2 = std::min(he2, ele.side(j)->measure()); break;
1047  case 3: he3 = std::min(he3, ele.side(j)->measure()); break;
1048  }
1049  }
1050  }
1051  if(h1 == d_max) h1 = 0;
1052  if(h2 == d_max) h2 = 0;
1053  if(h3 == d_max) h3 = 0;
1054  if(he2 == d_max) he2 = 0;
1055  if(he3 == d_max) he3 = 0;
1056 
1057  FILE * file;
1058  file = fopen(output_file.c_str(),"a");
1059  fprintf(file, "nA = %d;\n", mh_dh.side_ds->size());
1060  fprintf(file, "nB = %d;\n", mh_dh.el_ds->size());
1061  fprintf(file, "nBF = %d;\n", mh_dh.edge_ds->size());
1062  fprintf(file, "h1 = %e;\nh2 = %e;\nh3 = %e;\n", h1, h2, h3);
1063  fprintf(file, "he2 = %e;\nhe3 = %e;\n", he2, he3);
1064  fclose(file);
1065 }
1066 
1067 
1069  START_TIMER("DarcyFlowMH_Steady::set_mesh_data_for_bddc");
1070  // prepare mesh for BDDCML
1071  // initialize arrays
1072  // auxiliary map for creating coordinates of local dofs and global-to-local numbering
1073  std::map<int,arma::vec3> localDofMap;
1074  // connectivity for the subdomain, i.e. global dof numbers on element, stored element-by-element
1075  // Indices of Nodes on Elements
1076  std::vector<int> inet;
1077  // number of degrees of freedom on elements - determines elementwise chunks of INET array
1078  // Number of Nodes on Elements
1079  std::vector<int> nnet;
1080  // Indices of Subdomain Elements in Global Numbering - for local elements, their global indices
1081  std::vector<int> isegn;
1082 
1083  // This array is currently not used in BDDCML, it was used as an interface scaling alternative to scaling
1084  // by diagonal. It corresponds to the rho-scaling.
1085  std::vector<double> element_permeability;
1086 
1087  // maximal and minimal dimension of elements
1088  uint elDimMax = 1;
1089  uint elDimMin = 3;
1090  for ( unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++ ) {
1091  auto ele_ac = mh_dh.accessor(i_loc);
1092  // for each element, create local numbering of dofs as fluxes (sides), pressure (element centre), Lagrange multipliers (edges), compatible connections
1093 
1094  elDimMax = std::max( elDimMax, ele_ac.dim() );
1095  elDimMin = std::min( elDimMin, ele_ac.dim() );
1096 
1097  isegn.push_back( ele_ac.ele_global_idx() );
1098  int nne = 0;
1099 
1100  for (unsigned int si=0; si<ele_ac.element_accessor()->n_sides(); si++) {
1101  // insert local side dof
1102  int side_row = ele_ac.side_row(si);
1103  arma::vec3 coord = ele_ac.side(si)->centre();
1104 
1105  localDofMap.insert( std::make_pair( side_row, coord ) );
1106  inet.push_back( side_row );
1107  nne++;
1108  }
1109 
1110  // insert local pressure dof
1111  int el_row = ele_ac.ele_row();
1112  arma::vec3 coord = ele_ac.centre();
1113  localDofMap.insert( std::make_pair( el_row, coord ) );
1114  inet.push_back( el_row );
1115  nne++;
1116 
1117  for (unsigned int si=0; si<ele_ac.element_accessor()->n_sides(); si++) {
1118  // insert local edge dof
1119  int edge_row = ele_ac.edge_row(si);
1120  arma::vec3 coord = ele_ac.side(si)->centre();
1121 
1122  localDofMap.insert( std::make_pair( edge_row, coord ) );
1123  inet.push_back( edge_row );
1124  nne++;
1125  }
1126 
1127  // insert dofs related to compatible connections
1128  for ( unsigned int i_neigh = 0; i_neigh < ele_ac.element_accessor()->n_neighs_vb(); i_neigh++) {
1129  int edge_row = mh_dh.row_4_edge[ ele_ac.element_accessor()->neigh_vb[i_neigh]->edge_idx() ];
1130  arma::vec3 coord = ele_ac.element_accessor()->neigh_vb[i_neigh]->edge()->side(0)->centre();
1131 
1132  localDofMap.insert( std::make_pair( edge_row, coord ) );
1133  inet.push_back( edge_row );
1134  nne++;
1135  }
1136 
1137  nnet.push_back( nne );
1138 
1139  // version for rho scaling
1140  // trace computation
1141  arma::vec3 centre = ele_ac.centre();
1142  double conduct = data_->conductivity.value( centre , ele_ac.element_accessor() );
1143  auto aniso = data_->anisotropy.value( centre, ele_ac.element_accessor() );
1144 
1145  // compute mean on the diagonal
1146  double coef = 0.;
1147  for ( int i = 0; i < 3; i++) {
1148  coef = coef + aniso.at(i,i);
1149  }
1150  // Maybe divide by cs
1151  coef = conduct*coef / 3;
1152 
1153  OLD_ASSERT( coef > 0.,
1154  "Zero coefficient of hydrodynamic resistance %f . \n ", coef );
1155  element_permeability.push_back( 1. / coef );
1156  }
1157  //convert set of dofs to vectors
1158  // number of nodes (= dofs) on the subdomain
1159  int numNodeSub = localDofMap.size();
1160  OLD_ASSERT_EQUAL( (unsigned int)numNodeSub, mh_dh.global_row_4_sub_row->size() );
1161  // Indices of Subdomain Nodes in Global Numbering - for local nodes, their global indices
1162  std::vector<int> isngn( numNodeSub );
1163  // pseudo-coordinates of local nodes (i.e. dofs)
1164  // they need not be exact, they are used just for some geometrical considerations in BDDCML,
1165  // such as selection of corners maximizing area of a triangle, bounding boxes fro subdomains to
1166  // find candidate neighbours etc.
1167  std::vector<double> xyz( numNodeSub * 3 ) ;
1168  int ind = 0;
1169  std::map<int,arma::vec3>::iterator itB = localDofMap.begin();
1170  for ( ; itB != localDofMap.end(); ++itB ) {
1171  isngn[ind] = itB -> first;
1172 
1173  arma::vec3 coord = itB -> second;
1174  for ( int j = 0; j < 3; j++ ) {
1175  xyz[ j*numNodeSub + ind ] = coord[j];
1176  }
1177 
1178  ind++;
1179  }
1180  localDofMap.clear();
1181 
1182  // Number of Nodal Degrees of Freedom
1183  // nndf is trivially one - dofs coincide with nodes
1184  std::vector<int> nndf( numNodeSub, 1 );
1185 
1186  // prepare auxiliary map for renumbering nodes
1187  typedef std::map<int,int> Global2LocalMap_; //! type for storage of global to local map
1188  Global2LocalMap_ global2LocalNodeMap;
1189  for ( unsigned ind = 0; ind < isngn.size(); ++ind ) {
1190  global2LocalNodeMap.insert( std::make_pair( static_cast<unsigned>( isngn[ind] ), ind ) );
1191  }
1192 
1193  // renumber nodes in the inet array to locals
1194  int indInet = 0;
1195  for ( unsigned int iEle = 0; iEle < isegn.size(); iEle++ ) {
1196  int nne = nnet[ iEle ];
1197  for ( int ien = 0; ien < nne; ien++ ) {
1198 
1199  int indGlob = inet[indInet];
1200  // map it to local node
1201  Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1202  OLD_ASSERT( pos != global2LocalNodeMap.end(),
1203  "Cannot remap node index %d to local indices. \n ", indGlob );
1204  int indLoc = static_cast<int> ( pos -> second );
1205 
1206  // store the node
1207  inet[ indInet++ ] = indLoc;
1208  }
1209  }
1210 
1211  int numNodes = size;
1212  int numDofsInt = size;
1213  int spaceDim = 3; // TODO: what is the proper value here?
1214  int meshDim = elDimMax;
1215 
1216  bddc_ls -> load_mesh( spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1217 }
1218 
1219 
1220 
1221 
1222 //=============================================================================
1223 // DESTROY WATER MH SYSTEM STRUCTURE
1224 //=============================================================================
1226  if (previous_solution != nullptr) chkerr(VecDestroy(&previous_solution));
1227  if (steady_diagonal != nullptr) chkerr(VecDestroy(&steady_diagonal));
1228  if (new_diagonal != nullptr) chkerr(VecDestroy(&new_diagonal));
1229  if (steady_rhs != nullptr) chkerr(VecDestroy(&steady_rhs));
1230 
1231 
1232  if (schur0 != NULL) {
1233  delete schur0;
1234  chkerr(VecScatterDestroy(&par_to_all));
1235  }
1236 
1237  if (solution != NULL) {
1238  chkerr(VecDestroy(&sol_vec));
1239  delete [] solution;
1240  }
1241 
1242  if (output_object) delete output_object;
1243 
1244  if(time_ != nullptr)
1245  delete time_;
1246 
1247 }
1248 
1249 
1250 // ================================================
1251 // PARALLLEL PART
1252 //
1253 
1254 
1256  START_TIMER("prepare scatter");
1257  // prepare Scatter form parallel to sequantial in original numbering
1258  {
1259  IS is_loc;
1260  int i, *loc_idx; //, si;
1261 
1262  // create local solution vector
1263  solution = new double[size];
1264  VecCreateSeqWithArray(PETSC_COMM_SELF,1, size, solution,
1265  &(sol_vec));
1266 
1267  // create seq. IS to scatter par solutin to seq. vec. in original order
1268  // use essentialy row_4_id arrays
1269  loc_idx = new int [size];
1270  i = 0;
1271  for (auto ele : mesh_->elements_range()) {
1272  for (unsigned int si=0; si<ele->n_sides(); si++) {
1273  loc_idx[i++] = mh_dh.side_row_4_id[ mh_dh.side_dof( ele.side(si) ) ];
1274  }
1275  }
1276  for (auto ele : mesh_->elements_range()) {
1277  loc_idx[i++] = mh_dh.row_4_el[ ele.idx() ];
1278  }
1279  for(unsigned int i_edg=0; i_edg < mesh_->n_edges(); i_edg++) {
1280  loc_idx[i++] = mh_dh.row_4_edge[i_edg];
1281  }
1282  OLD_ASSERT( i==size,"Size of array does not match number of fills.\n");
1283  //DBGPRINT_INT("loc_idx",size,loc_idx);
1284  ISCreateGeneral(PETSC_COMM_SELF, size, loc_idx, PETSC_COPY_VALUES, &(is_loc));
1285  delete [] loc_idx;
1286  VecScatterCreate(schur0->get_solution(), is_loc, sol_vec,
1287  PETSC_NULL, &par_to_all);
1288  chkerr(ISDestroy(&(is_loc)));
1289  }
1291 
1292  END_TIMER("prepare scatter");
1293 
1294 }
1295 
1296 
1297 /*
1298 void mat_count_off_proc_values(Mat m, Vec v) {
1299  int n, first, last;
1300  const PetscInt *cols;
1301  Distribution distr(v);
1302 
1303  int n_off = 0;
1304  int n_on = 0;
1305  int n_off_rows = 0;
1306  MatGetOwnershipRange(m, &first, &last);
1307  for (int row = first; row < last; row++) {
1308  MatGetRow(m, row, &n, &cols, PETSC_NULL);
1309  bool exists_off = false;
1310  for (int i = 0; i < n; i++)
1311  if (distr.get_proc(cols[i]) != distr.myp())
1312  n_off++, exists_off = true;
1313  else
1314  n_on++;
1315  if (exists_off)
1316  n_off_rows++;
1317  MatRestoreRow(m, row, &n, &cols, PETSC_NULL);
1318  }
1319 }
1320 */
1321 
1322 
1323 
1324 
1325 
1326 
1327 
1328 
1329 
1330 
1331 
1332 
1334 {
1335  double *local_sol = schur0->get_solution_array();
1336 
1337  // cycle over local element rows
1338 
1339  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1340  for (unsigned int i_loc_el = 0; i_loc_el < mh_dh.el_ds->lsize(); i_loc_el++) {
1341  auto ele_ac = mh_dh.accessor(i_loc_el);
1342  // set initial condition
1343  local_sol[ele_ac.ele_local_row()] = data_->init_pressure.value(ele_ac.centre(),ele_ac.element_accessor());
1344  }
1345 
1347 
1348 }
1349 
1351  // save diagonal of steady matrix
1352  MatGetDiagonal(*( schur0->get_matrix() ), steady_diagonal);
1353  // save RHS
1354  VecCopy(*( schur0->get_rhs()), steady_rhs);
1355 
1356 
1357  PetscScalar *local_diagonal;
1358  VecGetArray(new_diagonal,& local_diagonal);
1359 
1360  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1361 
1362  balance_->start_mass_assembly(data_->water_balance_idx);
1363 
1364  //DebugOut().fmt("time_term lsize: {} {}\n", mh_dh.el_ds->myp(), mh_dh.el_ds->lsize());
1365  for (unsigned int i_loc_el = 0; i_loc_el < mh_dh.el_ds->lsize(); i_loc_el++) {
1366  auto ele_ac = mh_dh.accessor(i_loc_el);
1367 
1368  // set new diagonal
1369  double diagonal_coeff = data_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor())
1370  * data_->storativity.value(ele_ac.centre(), ele_ac.element_accessor())
1371  * ele_ac.measure();
1372  local_diagonal[ele_ac.ele_local_row()]= - diagonal_coeff / time_->dt();
1373 
1374  //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);
1375  balance_->add_mass_matrix_values(data_->water_balance_idx,
1376  ele_ac.region().bulk_idx(), { LongIdx(ele_ac.ele_row()) }, {diagonal_coeff});
1377  }
1378  VecRestoreArray(new_diagonal,& local_diagonal);
1379  MatDiagonalSet(*( schur0->get_matrix() ), new_diagonal, ADD_VALUES);
1380 
1383 
1384  balance_->finish_mass_assembly(data_->water_balance_idx);
1385 }
1386 
1388  START_TIMER("modify system");
1389  if (time_->is_changed_dt() && time_->step().index()>0) {
1390  double scale_factor=time_->step(-2).length()/time_->step().length();
1391  if (scale_factor != 1.0) {
1392  // if time step has changed and setup_time_term not called
1393  MatDiagonalSet(*( schur0->get_matrix() ),steady_diagonal, INSERT_VALUES);
1394 
1395  VecScale(new_diagonal, time_->last_dt()/time_->dt());
1396  MatDiagonalSet(*( schur0->get_matrix() ),new_diagonal, ADD_VALUES);
1398  }
1399  }
1400 
1401  // modify RHS - add previous solution
1402  //VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, previous_solution);
1403  VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, schur0->get_solution());
1404  VecAXPY(*( schur0->get_rhs()), 1.0, steady_rhs);
1406 
1407  //VecSwap(previous_solution, schur0->get_solution());
1408 }
1409 
1410 
1411 //-----------------------------------------------------------------------------
1412 // vim: set cindent:
TimeGovernor & time()
Definition: equation.hh:148
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
unsigned int size() const
get global size
void get_solution_vector(double *&vec, unsigned int &vec_size) override
void make_serial_scatter()
FieldSet * eq_data_
Definition: equation.hh:232
Distribution * side_ds
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:49
void reinit(Mesh *mesh)
Output class for darcy_flow_mh model.
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:329
virtual void initialize_specific()
void get_parallel_solution_vector(Vec &vector) override
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
SchurComplement SchurComplement
Definition: linsys.hh:92
void assembly_mh_matrix(MultidimAssembly &assembler)
void set_symmetric(bool flag=true)
Definition: linsys.hh:541
unsigned int uint
Classes with algorithms for computation of intersections of meshes.
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void init_eq_data()
virtual ~DarcyMH() override
MixedMeshIntersections & mixed_intersections()
Definition: mesh.cc:710
static const Input::Type::Record & get_input_type()
The specification of output stream.
Definition: output_time.cc:38
Common base for intersection object.
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Class for declaration of the input of type Bool.
Definition: type_base.hh:459
void set_from_input(const Input::Record in_rec) override
unsigned int edge_idx() const
Definition: side_impl.hh:62
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:243
virtual void start_add_assembly()
Definition: linsys.hh:341
virtual void output_data() override
Write computed fields.
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:264
Wrappers for linear systems based on MPIAIJ and MATIS format.
bool is_end() const
Returns true if the actual time is greater than or equal to the end time.
void next_time()
Proceed to the next time according to current estimated time step.
void set_rhs_changed()
Definition: linsys.hh:218
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:110
Vec steady_rhs
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
void initialize() override
LongIdx * side_row_4_id
Definition: mesh.h:76
static const Input::Type::Record & type_field_descriptor()
virtual void start_allocation()
Definition: linsys.hh:333
Distribution * el_ds
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:147
Mat< double, N, 1 > vec
Definition: armor.hh:211
virtual void finish_assembly()=0
virtual double get_solution_precision()=0
std::vector< std::vector< ILpair > > element_intersections_
SideIter side(const unsigned int loc_index)
Definition: accessors.hh:137
const RegionDB & region_db() const
Definition: mesh.h:143
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
const TimeStep & step(int index=-1) const
static const std::string field_descriptor_record_description(const string &record_name)
Definition: field_common.cc:73
Class for declaration of the integral input data.
Definition: type_base.hh:490
Basic time management functionality for unsteady (and steady) solvers (class Equation).
unsigned int bulk_ele_idx() const
Returns index of bulk element.
LongIdx * row_4_el
bool solution_changed_for_scatter
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
static const Input::Type::Record & get_input_type()
Definition: linsys_BDDC.cc:37
std::shared_ptr< EqData > data_
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
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:346
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:727
void view(const char *name="") const
std::shared_ptr< LocalToGlobalMap > global_row_4_sub_row
Necessary only for BDDC solver.
unsigned int min_n_it_
void modify_system()
unsigned int nonlinear_iteration_
void set_mesh_data_for_bddc(LinSys_BDDC *bddc_ls)
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
void update_solution() override
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void set_from_input(const Input::Record in_rec) override
Definition: linsys_BDDC.cc:302
virtual void postprocess()
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:541
static const Input::Type::Instance & get_input_type_specific()
double * get_solution_array()
Definition: linsys.hh:325
std::shared_ptr< Balance > balance_
#define MPI_MIN
Definition: mpi.h:198
unsigned int edge_idx()
Definition: neighbours.h:153
FLOW123D_FORCE_LINK_IN_CHILD(darcy_flow_mh)
int n_schur_compls
const Vec & get_solution()
Definition: linsys.hh:282
Input::Type::Record type() const
Definition: accessors.cc:273
FMT_FUNC int fprintf(std::ostream &os, CStringRef format, ArgList args)
Definition: ostream.cc:56
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Ret val(const string &key) const
unsigned int n_sides() const
Definition: elements.h:135
#define xprintf(...)
Definition: system.hh:92
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
Distribution * edge_ds
Mesh * mesh_
Definition: equation.hh:223
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
double length() const
unsigned int side_dof(const SideIter side) const
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:503
static Input::Type::Abstract & get_input_type()
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1041
bool use_steady_assembly_
unsigned int n_sides() const
Definition: mesh.cc:227
void allocate_mh_matrix()
static const Input::Type::Selection & get_bc_type_selection()
Return a Selection corresponding to enum BC_Type.
virtual void assembly_linear_system()
virtual const Vec * get_rhs()
Definition: linsys.hh:203
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
LinSys * schur0
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:273
void set_solution(Vec sol_vec)
Definition: linsys.hh:290
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
DarcyMH(Mesh &mesh, const Input::Record in_rec)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
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 create_linear_system(Input::AbstractRecord rec)
Initialize global_row_4_sub_row.
MH_DofHandler mh_dh
LongIdx * row_4_edge
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:346
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.
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
static const Input::Type::Record & get_input_type()
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
virtual SolveInfo solve()=0
#define END_TIMER(tag)
Ends a timer with specified tag.
#define OLD_ASSERT_EQUAL(a, b)
Definition: global_defs.h:133
EqData()
Creation of all fields.
Definition: system.hh:64
void set_matrix_changed()
Definition: linsys.hh:212
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:32
void print_matlab_matrix(string matlab_file)
Print darcy flow matrix in matlab format into a file.
mixed-hybrid model of linear Darcy flow, possibly unsteady.
unsigned int index() const
void zero_time_step() override
void set_positive_definite(bool flag=true)
Definition: linsys.hh:556
virtual void read_initial_condition()
Record type proxy class.
Definition: type_record.hh:182
void prepare_parallel_bddc()
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:137
virtual void prepare_new_time_step()
postprocess velocity field (add sources)
Class for representation SI units of Fields.
Definition: unit_si.hh:40
virtual const Mat * get_matrix()
Definition: linsys.hh:187
double last_dt() const
Vec new_diagonal
virtual double solution_precision() const
virtual void setup_time_term()
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
#define DebugOut()
Macro defining &#39;debug&#39; record of log.
Definition: logger.hh:252
virtual double compute_residual()=0
Vec previous_solution
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
Template for classes storing finite set of named values.
Implementation of range helper class.
void rhs_set_value(int row, double val)
Definition: linsys.hh:381
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