Flow123d  JS_before_hm-1016-g09ec2cb
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 "system/index_types.hh"
35 #include "input/factory.hh"
36 
37 #include "mesh/mesh.h"
38 #include "mesh/partitioning.hh"
39 #include "mesh/accessors.hh"
40 #include "mesh/range_wrapper.hh"
41 #include "la/distribution.hh"
42 #include "la/linsys.hh"
43 #include "la/linsys_PETSC.hh"
44 #include "la/linsys_BDDC.hh"
45 #include "la/schur.hh"
46 //#include "la/sparse_graph.hh"
48 #include "la/vector_mpi.hh"
49 
50 #include "flow/assembly_mh.hh"
51 #include "flow/darcy_flow_mh.hh"
53 
54 #include "tools/time_governor.hh"
56 #include "fields/field.hh"
57 #include "fields/field_values.hh"
59 #include "fields/field_fe.hh"
60 #include "fields/field_divide.hh"
61 
62 #include "coupling/balance.hh"
63 
66 
67 #include "fem/fe_p.hh"
68 
69 
70 FLOW123D_FORCE_LINK_IN_CHILD(darcy_flow_mh)
71 
72 
73 
74 
75 namespace it = Input::Type;
76 
77 const it::Selection & DarcyMH::get_mh_mortar_selection() {
78  return it::Selection("MH_MortarMethod")
79  .add_value(NoMortar, "None", "No Mortar method is applied.")
80  .add_value(MortarP0, "P0", "Mortar space: P0 on elements of lower dimension.")
81  .add_value(MortarP1, "P1", "Mortar space: P1 on intersections, using non-conforming pressures.")
82  .close();
83 }
84 
85 
87  return it::Selection("Flow_Darcy_BC_Type")
88  .add_value(none, "none",
89  "Homogeneous Neumann boundary condition\n(zero normal flux over the boundary).")
90  .add_value(dirichlet, "dirichlet",
91  "Dirichlet boundary condition. "
92  "Specify the pressure head through the ``bc_pressure`` field "
93  "or the piezometric head through the ``bc_piezo_head`` field.")
94  .add_value(total_flux, "total_flux", "Flux boundary condition (combines Neumann and Robin type). "
95  "Water inflow equal to (($ \\delta_d(q_d^N + \\sigma_d (h_d^R - h_d) )$)). "
96  "Specify the water inflow by the ``bc_flux`` field, the transition coefficient by ``bc_robin_sigma`` "
97  "and the reference pressure head or piezometric head through ``bc_pressure`` or ``bc_piezo_head`` respectively.")
98  .add_value(seepage, "seepage",
99  "Seepage face boundary condition. Pressure and inflow bounded from above. Boundary with potential seepage flow "
100  "is described by the pair of inequalities: "
101  "(($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. "
102  "Caution: setting (($q_d^N$)) strictly negative "
103  "may lead to an ill posed problem since a positive outflow is enforced. "
104  "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."
105  )
106  .add_value(river, "river",
107  "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: "
108  "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d) )$)). For the water level under the bedrock, constant infiltration is used: "
109  "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d^S) )$)). Parameters: ``bc_pressure``, ``bc_switch_pressure``, "
110  " ``bc_sigma``, ``bc_flux``."
111  )
112  .close();
113 }
114 
116 
117  const it::Record &field_descriptor =
118  it::Record("Flow_Darcy_MH_Data",FieldCommon::field_descriptor_record_description("Flow_Darcy_MH_Data") )
119  .copy_keys( DarcyMH::EqData().make_field_descriptor_type("Flow_Darcy_MH_Data_aux") )
120  .declare_key("bc_piezo_head", FieldAlgorithmBase< 3, FieldValue<3>::Scalar >::get_input_type_instance(),
121  "Boundary piezometric head for BC types: dirichlet, robin, and river." )
122  .declare_key("bc_switch_piezo_head", FieldAlgorithmBase< 3, FieldValue<3>::Scalar >::get_input_type_instance(),
123  "Boundary switch piezometric head for BC types: seepage, river." )
124  .declare_key("init_piezo_head", FieldAlgorithmBase< 3, FieldValue<3>::Scalar >::get_input_type_instance(),
125  "Initial condition for the pressure given as the piezometric head." )
126  .close();
127  return field_descriptor;
128 }
129 
131 
132  it::Record ns_rec = Input::Type::Record("NonlinearSolver", "Non-linear solver settings.")
133  .declare_key("linear_solver", LinSys::get_input_type(), it::Default("{}"),
134  "Linear solver for MH problem.")
135  .declare_key("tolerance", it::Double(0.0), it::Default("1E-6"),
136  "Residual tolerance.")
137  .declare_key("min_it", it::Integer(0), it::Default("1"),
138  "Minimum number of iterations (linear solutions) to use.\nThis is usefull if the convergence criteria "
139  "does not characterize your goal well enough so it converges prematurely, possibly even without a single linear solution."
140  "If greater then 'max_it' the value is set to 'max_it'.")
141  .declare_key("max_it", it::Integer(0), it::Default("100"),
142  "Maximum number of iterations (linear solutions) of the non-linear solver.")
143  .declare_key("converge_on_stagnation", it::Bool(), it::Default("false"),
144  "If a stagnation of the nonlinear solver is detected the solver stops. "
145  "A divergence is reported by default, forcing the end of the simulation. By setting this flag to 'true', the solver "
146  "ends with convergence success on stagnation, but it reports warning about it.")
147  .close();
148 
149  DarcyMH::EqData eq_data;
150 
151  return it::Record("Flow_Darcy_MH", "Mixed-Hybrid solver for saturated Darcy flow.")
154  .declare_key("gravity", it::Array(it::Double(), 3,3), it::Default("[ 0, 0, -1]"),
155  "Vector of the gravity force. Dimensionless.")
157  "Input data for Darcy flow model.")
158  .declare_key("nonlinear_solver", ns_rec, it::Default("{}"),
159  "Non-linear solver for MH problem.")
160  .declare_key("output_stream", OutputTime::get_input_type(), it::Default("{}"),
161  "Output stream settings.\n Specify file format, precision etc.")
162 
163  .declare_key("output", DarcyFlowMHOutput::get_input_type(eq_data, "Flow_Darcy_MH"),
164  IT::Default("{ \"fields\": [ \"pressure_p0\", \"velocity_p0\" ] }"),
165  "Specification of output fields and output times.")
167  "Output settings specific to Darcy flow model.\n"
168  "Includes raw output and some experimental functionality.")
169  .declare_key("balance", Balance::get_input_type(), it::Default("{}"),
170  "Settings for computing mass balance.")
171  .declare_key("n_schurs", it::Integer(0,2), it::Default("2"),
172  "Number of Schur complements to perform when solving MH system.")
173  .declare_key("mortar_method", get_mh_mortar_selection(), it::Default("\"None\""),
174  "Method for coupling Darcy flow between dimensions on incompatible meshes. [Experimental]" )
175  .close();
176 }
177 
178 
179 const int DarcyMH::registrar =
180  Input::register_class< DarcyMH, Mesh &, const Input::Record >("Flow_Darcy_MH") +
182 
183 
184 
186 {
188 
189  *this += field_ele_pressure.name("pressure_p0")
190  .units(UnitSI().m())
192  .description("Pressure solution - P0 interpolation.");
193 
194  *this += field_edge_pressure.name("pressure_edge")
195  .units(UnitSI().m())
197  .description("Pressure solution - Crouzeix-Raviart interpolation.");
198 
199  *this += field_ele_piezo_head.name("piezo_head_p0")
200  .units(UnitSI().m())
202  .description("Piezo head solution - P0 interpolation.");
203 
204  *this += field_ele_velocity.name("velocity_p0")
205  .units(UnitSI().m().s(-1))
207  .description("Velocity solution - P0 interpolation.");
208 
209  *this += flux.name("flux")
210  .units(UnitSI().m().s(-1))
212  .description("Darcy flow flux.");
213 
214  *this += anisotropy.name("anisotropy")
215  .description("Anisotropy of the conductivity tensor.")
216  .input_default("1.0")
218 
219  *this += cross_section.name("cross_section")
220  .description("Complement dimension parameter (cross section for 1D, thickness for 2D).")
221  .input_default("1.0")
222  .units( UnitSI().m(3).md() );
223 
224  *this += conductivity.name("conductivity")
225  .description("Isotropic conductivity scalar.")
226  .input_default("1.0")
227  .units( UnitSI().m().s(-1) )
228  .set_limits(0.0);
229 
230  *this += sigma.name("sigma")
231  .description("Transition coefficient between dimensions.")
232  .input_default("1.0")
234 
235  *this += water_source_density.name("water_source_density")
236  .description("Water source density.")
237  .input_default("0.0")
238  .units( UnitSI().s(-1) );
239 
240  *this += bc_type.name("bc_type")
241  .description("Boundary condition type.")
243  .input_default("\"none\"")
245 
246  *this += bc_pressure
248  .name("bc_pressure")
249  .description("Prescribed pressure value on the boundary. Used for all values of ``bc_type`` except ``none`` and ``seepage``. "
250  "See documentation of ``bc_type`` for exact meaning of ``bc_pressure`` in individual boundary condition types.")
251  .input_default("0.0")
252  .units( UnitSI().m() );
253 
254  *this += bc_flux
256  .name("bc_flux")
257  .description("Incoming water boundary flux. Used for bc_types : ``total_flux``, ``seepage``, ``river``.")
258  .input_default("0.0")
259  .units( UnitSI().m().s(-1) );
260 
261  *this += bc_robin_sigma
263  .name("bc_robin_sigma")
264  .description("Conductivity coefficient in the ``total_flux`` or the ``river`` boundary condition type.")
265  .input_default("0.0")
266  .units( UnitSI().s(-1) );
267 
268  *this += bc_switch_pressure
270  .name("bc_switch_pressure")
271  .description("Critical switch pressure for ``seepage`` and ``river`` boundary conditions.")
272  .input_default("0.0")
273  .units( UnitSI().m() );
274 
275 
276  //these are for unsteady
277  *this += init_pressure.name("init_pressure")
278  .description("Initial condition for pressure in time dependent problems.")
279  .input_default("0.0")
280  .units( UnitSI().m() );
281 
282  *this += storativity.name("storativity")
283  .description("Storativity (in time dependent problems).")
284  .input_default("0.0")
285  .units( UnitSI().m(-1) );
286 
287  *this += extra_storativity.name("extra_storativity")
288  .description("Storativity added from upstream equation.")
289  .units( UnitSI().m(-1) )
290  .input_default("0.0")
291  .flags( input_copy );
292 
293  *this += extra_source.name("extra_water_source_density")
294  .description("Water source density added from upstream equation.")
295  .input_default("0.0")
296  .units( UnitSI().s(-1) )
297  .flags( input_copy );
298 
299  //time_term_fields = this->subset({"storativity"});
300  //main_matrix_fields = this->subset({"anisotropy", "conductivity", "cross_section", "sigma", "bc_type", "bc_robin_sigma"});
301  //rhs_fields = this->subset({"water_source_density", "bc_pressure", "bc_flux"});
302 }
303 
304 
305 
306 //=============================================================================
307 // CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
308 // - do it in parallel:
309 // - initial distribution of elements, edges
310 //
311 /*! @brief CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
312  *
313  * Parameters {Solver,NSchurs} number of performed Schur
314  * complements (0,1,2) for water flow MH-system
315  *
316  */
317 //=============================================================================
318 DarcyMH::DarcyMH(Mesh &mesh_in, const Input::Record in_rec, TimeGovernor *tm)
319 : DarcyFlowInterface(mesh_in, in_rec),
320  output_object(nullptr),
321  data_changed_(false),
322  schur0(nullptr),
323  steady_diagonal(nullptr),
324  steady_rhs(nullptr),
325  new_diagonal(nullptr),
326  previous_solution(nullptr)
327 {
328 
329  START_TIMER("Darcy constructor");
330  {
331  auto time_rec = in_rec.val<Input::Record>("time");
332  if (tm == nullptr)
333  {
334  time_ = new TimeGovernor(time_rec);
335  }
336  else
337  {
338  TimeGovernor tm_from_rec(time_rec);
339  if (!tm_from_rec.is_default()) // is_default() == false when time record is present in input file
340  {
341  MessageOut() << "Duplicate key 'time', time in flow equation is already initialized from parent class!";
342  ASSERT(false);
343  }
344  time_ = tm;
345  }
346  }
347 
348  data_ = make_shared<EqData>();
350 
351  data_->is_linear=true;
352 
353  size = mesh_->n_elements() + mesh_->n_sides() + mesh_->n_edges();
354  n_schur_compls = in_rec.val<int>("n_schurs");
355  data_->mortar_method_= in_rec.val<MortarMethod>("mortar_method");
356  if (data_->mortar_method_ != NoMortar) {
358  }
359 
360 
361 
362  //side_ds->view( std::cout );
363  //el_ds->view( std::cout );
364  //edge_ds->view( std::cout );
365  //rows_ds->view( std::cout );
366 
367 }
368 
369 
370 
372 //connecting data fields with mesh
373 {
374 
375  START_TIMER("data init");
376  data_->mesh = mesh_;
377  data_->set_mesh(*mesh_);
378 
379  auto gravity_array = input_record_.val<Input::Array>("gravity");
380  std::vector<double> gvec;
381  gravity_array.copy_to(gvec);
382  gvec.push_back(0.0); // zero pressure shift
383  data_->gravity_ = arma::vec(gvec);
384  data_->gravity_vec_ = data_->gravity_.subvec(0,2);
385 
386  data_->bc_pressure.add_factory(
387  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
388  (data_->gravity_, "bc_piezo_head") );
389  data_->bc_switch_pressure.add_factory(
390  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
391  (data_->gravity_, "bc_switch_piezo_head") );
392  data_->init_pressure.add_factory(
393  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
394  (data_->gravity_, "init_piezo_head") );
395 
396 
397  data_->set_input_list( this->input_record_.val<Input::Array>("input_fields"), *time_ );
398  // Check that the time step was set for the transient simulation.
399  if (! zero_time_term(true) && time_->is_default() ) {
400  //THROW(ExcAssertMsg());
401  //THROW(ExcMissingTimeGovernor() << input_record_.ei_address());
402  MessageOut() << "Missing the key 'time', obligatory for the transient problems." << endl;
403  ASSERT(false);
404  }
405 
406  data_->mark_input_times(*time_);
407 }
408 
409 
410 
412 
413  { // init DOF handler for pressure fields
414 // std::shared_ptr< FiniteElement<0> > fe0_rt = std::make_shared<FE_RT0_disc<0>>();
415  std::shared_ptr< FiniteElement<1> > fe1_rt = std::make_shared<FE_RT0_disc<1>>();
416  std::shared_ptr< FiniteElement<2> > fe2_rt = std::make_shared<FE_RT0_disc<2>>();
417  std::shared_ptr< FiniteElement<3> > fe3_rt = std::make_shared<FE_RT0_disc<3>>();
418  std::shared_ptr< FiniteElement<0> > fe0_disc = std::make_shared<FE_P_disc<0>>(0);
419  std::shared_ptr< FiniteElement<1> > fe1_disc = std::make_shared<FE_P_disc<1>>(0);
420  std::shared_ptr< FiniteElement<2> > fe2_disc = std::make_shared<FE_P_disc<2>>(0);
421  std::shared_ptr< FiniteElement<3> > fe3_disc = std::make_shared<FE_P_disc<3>>(0);
422  std::shared_ptr< FiniteElement<0> > fe0_cr = std::make_shared<FE_CR<0>>();
423  std::shared_ptr< FiniteElement<1> > fe1_cr = std::make_shared<FE_CR<1>>();
424  std::shared_ptr< FiniteElement<2> > fe2_cr = std::make_shared<FE_CR<2>>();
425  std::shared_ptr< FiniteElement<3> > fe3_cr = std::make_shared<FE_CR<3>>();
426 // static FiniteElement<0> fe0_sys = FE_P_disc<0>(0); //TODO fix and use solution with FESystem<0>( {fe0_rt, fe0_disc, fe0_cr} )
427  FESystem<0> fe0_sys( {fe0_disc, fe0_disc, fe0_cr} );
428  FESystem<1> fe1_sys( {fe1_rt, fe1_disc, fe1_cr} );
429  FESystem<2> fe2_sys( {fe2_rt, fe2_disc, fe2_cr} );
430  FESystem<3> fe3_sys( {fe3_rt, fe3_disc, fe3_cr} );
431  MixedPtr<FESystem> fe_sys( std::make_shared<FESystem<0>>(fe0_sys), std::make_shared<FESystem<1>>(fe1_sys),
432  std::make_shared<FESystem<2>>(fe2_sys), std::make_shared<FESystem<3>>(fe3_sys) );
433  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_sys);
434  data_->dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
435  data_->dh_->distribute_dofs(ds);
436  }
437 
438  init_eq_data();
439  this->data_->multidim_assembler = AssemblyBase::create< AssemblyMH >(data_);
441 
442  { // construct pressure, velocity and piezo head fields
443  uint rt_component = 0;
444  auto ele_flux_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(data_->dh_, rt_component);
445  data_->full_solution = ele_flux_ptr->vec();
446  data_->flux.set_field(mesh_->region_db().get_region_set("ALL"), ele_flux_ptr);
447 
448  auto ele_velocity_ptr = std::make_shared< FieldDivide<3, FieldValue<3>::VectorFixed> >(ele_flux_ptr, data_->cross_section);
449  data_->field_ele_velocity.set_field(mesh_->region_db().get_region_set("ALL"), ele_velocity_ptr);
450 
451  uint p_ele_component = 0;
452  auto ele_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(data_->dh_, p_ele_component, &data_->full_solution);
453  data_->field_ele_pressure.set_field(mesh_->region_db().get_region_set("ALL"), ele_pressure_ptr);
454 
455  uint p_edge_component = 1;
456  auto edge_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(data_->dh_, p_edge_component, &data_->full_solution);
457  data_->field_edge_pressure.set_field(mesh_->region_db().get_region_set("ALL"), edge_pressure_ptr);
458 
459  arma::vec4 gravity = (-1) * data_->gravity_;
460  auto ele_piezo_head_ptr = std::make_shared< FieldAddPotential<3, FieldValue<3>::Scalar> >(gravity, ele_pressure_ptr);
461  data_->field_ele_piezo_head.set_field(mesh_->region_db().get_region_set("ALL"), ele_piezo_head_ptr);
462  }
463 
464  { // init DOF handlers represents edge DOFs
465  uint p_edge_component = 2;
466  data_->dh_cr_ = std::make_shared<SubDOFHandlerMultiDim>(data_->dh_,p_edge_component);
467  }
468 
469  { // init DOF handlers represents side DOFs
470  MixedPtr<FE_CR_disc> fe_cr_disc;
471  std::shared_ptr<DiscreteSpace> ds_cr_disc = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_cr_disc);
472  data_->dh_cr_disc_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
473  data_->dh_cr_disc_->distribute_dofs(ds_cr_disc);
474  }
475 
476  // Initialize bc_switch_dirichlet to size of global boundary.
477  data_->bc_switch_dirichlet.resize(mesh_->n_elements()+mesh_->n_elements(true), 1);
478 
479 
482  .val<Input::Record>("nonlinear_solver")
483  .val<Input::AbstractRecord>("linear_solver");
484 
485  // auxiliary set_time call since allocation assembly evaluates fields as well
488 
489 
490 
491  // allocate time term vectors
492  VecDuplicate(schur0->get_solution(), &previous_solution);
493  VecCreateMPI(PETSC_COMM_WORLD, data_->dh_->distr()->lsize(),PETSC_DETERMINE,&(steady_diagonal));
494  VecDuplicate(steady_diagonal,& new_diagonal);
495  VecZeroEntries(new_diagonal);
496  VecDuplicate(steady_diagonal, &steady_rhs);
497 
498 
499  // initialization of balance object
500  balance_ = std::make_shared<Balance>("water", mesh_);
501  balance_->init_from_input(input_record_.val<Input::Record>("balance"), time());
502  data_->water_balance_idx = balance_->add_quantity("water_volume");
503  balance_->allocate(data_->dh_->distr()->lsize(), 1);
504  balance_->units(UnitSI().m(3));
505 
506 
507  data_->balance = balance_;
508  data_->lin_sys = schur0;
509 
510 
512 }
513 
515 {}
516 
518 {
519 
520  /* TODO:
521  * - Allow solution reconstruction (pressure and velocity) from initial condition on user request.
522  * - Steady solution as an intitial condition may be forced by setting inti_time =-1, and set data for the steady solver in that time.
523  * Solver should be able to switch from and to steady case depending on the zero time term.
524  */
525 
527 
528  // zero_time_term means steady case
529  bool zero_time_term_from_right = zero_time_term();
530 
531 
532  if (zero_time_term_from_right) {
533  // steady case
534  VecZeroEntries(schur0->get_solution());
535  //read_initial_condition(); // Possible solution guess for steady case.
536  use_steady_assembly_ = true;
537  solve_nonlinear(); // with right limit data
538  } else {
539  VecZeroEntries(schur0->get_solution());
540  VecZeroEntries(previous_solution);
542  assembly_linear_system(); // in particular due to balance
543 // print_matlab_matrix("matrix_zero");
544  // TODO: reconstruction of solution in zero time.
545  }
546  //solution_output(T,right_limit); // data for time T in any case
547  output_data();
548 }
549 
550 //=============================================================================
551 // COMPOSE and SOLVE WATER MH System possibly through Schur complements
552 //=============================================================================
554 {
555  START_TIMER("Solving MH system");
556 
557  time_->next_time();
558 
559  time_->view("DARCY"); //time governor information output
560 
561  solve_time_step();
562 
563  data_->full_solution.local_to_ghost_begin();
564  data_->full_solution.local_to_ghost_end();
565 }
566 
567 
569 {
571 
572  bool zero_time_term_from_left=zero_time_term();
573 
574  bool jump_time = data_->storativity.is_jump_time();
575  if (! zero_time_term_from_left) {
576  // time term not treated as zero
577  // Unsteady solution up to the T.
578 
579  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
580  use_steady_assembly_ = false;
581 
583  solve_nonlinear(); // with left limit data
584  if (jump_time) {
585  WarningOut() << "Output of solution discontinuous in time not supported yet.\n";
586  //solution_output(T, left_limit); // output use time T- delta*dt
587  //output_data();
588  }
589  }
590 
591  if (time_->is_end()) {
592  // output for unsteady case, end_time should not be the jump time
593  // but rether check that
594  if (! zero_time_term_from_left && ! jump_time && output) output_data();
595  return;
596  }
597 
599  bool zero_time_term_from_right=zero_time_term();
600  if (zero_time_term_from_right) {
601  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
602  use_steady_assembly_ = true;
603  solve_nonlinear(); // with right limit data
604 
605  } else if (! zero_time_term_from_left && jump_time) {
606  WarningOut() << "Discontinuous time term not supported yet.\n";
607  //solution_transfer(); // internally call set_time(T, left) and set_time(T,right) again
608  //solve_nonlinear(); // with right limit data
609  }
610 
611  //solution_output(T,right_limit); // data for time T in any case
612  if (output) output_data();
613 }
614 
615 
616 bool DarcyMH::zero_time_term(bool time_global) {
617  if (time_global) {
618  return (data_->storativity.input_list_size() == 0);
619  } else {
620  return data_->storativity.field_result(mesh_->region_db().get_region_set("BULK")) == result_zeros;
621  }
622 }
623 
624 
626 {
627 
629  double residual_norm = schur0->compute_residual();
631  MessageOut().fmt("[nonlinear solver] norm of initial residual: {}\n", residual_norm);
632 
633  // Reduce is_linear flag.
634  int is_linear_common;
635  MPI_Allreduce(&(data_->is_linear), &is_linear_common,1, MPI_INT ,MPI_MIN,PETSC_COMM_WORLD);
636 
637  Input::Record nl_solver_rec = input_record_.val<Input::Record>("nonlinear_solver");
638  this->tolerance_ = nl_solver_rec.val<double>("tolerance");
639  this->max_n_it_ = nl_solver_rec.val<unsigned int>("max_it");
640  this->min_n_it_ = nl_solver_rec.val<unsigned int>("min_it");
641  if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->max_n_it_;
642 
643  if (! is_linear_common) {
644  // set tolerances of the linear solver unless they are set by user.
645  schur0->set_tolerances(0.1*this->tolerance_, 0.01*this->tolerance_, 100);
646  }
647  vector<double> convergence_history;
648 
649  Vec save_solution;
650  VecDuplicate(schur0->get_solution(), &save_solution);
651  while (nonlinear_iteration_ < this->min_n_it_ ||
652  (residual_norm > this->tolerance_ && nonlinear_iteration_ < this->max_n_it_ )) {
653  OLD_ASSERT_EQUAL( convergence_history.size(), nonlinear_iteration_ );
654  convergence_history.push_back(residual_norm);
655 
656  // stagnation test
657  if (convergence_history.size() >= 5 &&
658  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
659  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
660  // stagnation
661  if (input_record_.val<Input::Record>("nonlinear_solver").val<bool>("converge_on_stagnation")) {
662  WarningOut().fmt("Accept solution on stagnation. Its: {} Residual: {}\n", nonlinear_iteration_, residual_norm);
663  break;
664  } else {
665  THROW(ExcSolverDiverge() << EI_Reason("Stagnation."));
666  }
667  }
668 
669  if (! is_linear_common)
670  VecCopy( schur0->get_solution(), save_solution);
673 
674  // hack to make BDDC work with empty compute_residual
675  if (is_linear_common){
676  // we want to print this info in linear (and steady) case
677  residual_norm = schur0->compute_residual();
678  MessageOut().fmt("[nonlinear solver] lin. it: {}, reason: {}, residual: {}\n",
679  si.n_iterations, si.converged_reason, residual_norm);
680  break;
681  }
682  data_changed_=true; // force reassembly for non-linear case
683 
684  double alpha = 1; // how much of new solution
685  VecAXPBY(schur0->get_solution(), (1-alpha), alpha, save_solution);
686 
687  //LogOut().fmt("Linear solver ended with reason: {} \n", si.converged_reason );
688  //OLD_ASSERT( si.converged_reason >= 0, "Linear solver failed to converge. Convergence reason %d \n", si.converged_reason );
690 
691  residual_norm = schur0->compute_residual();
692  MessageOut().fmt("[nonlinear solver] it: {} lin. it: {}, reason: {}, residual: {}\n",
693  nonlinear_iteration_, si.n_iterations, si.converged_reason, residual_norm);
694  }
695  chkerr(VecDestroy(&save_solution));
696  this -> postprocess();
697 
698  // adapt timestep
699  if (! this->zero_time_term()) {
700  double mult = 1.0;
701  if (nonlinear_iteration_ < 3) mult = 1.6;
702  if (nonlinear_iteration_ > 7) mult = 0.7;
703  time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
704  // int result = time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
705  //DebugOut().fmt("time adaptivity, res: {} it: {} m: {} dt: {} edt: {}\n", result, nonlinear_iteration_, mult, time_->dt(), time_->estimate_dt());
706  }
707 }
708 
710 {
712 }
713 
715 {
716  START_TIMER("postprocess");
717 
718  //fix velocity when mortar method is used
719  if(data_->mortar_method_ != MortarMethod::NoMortar){
720  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(data_);
721  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
722  unsigned int dim = dh_cell.dim();
723  multidim_assembler[dim-1]->fix_velocity(dh_cell);
724  }
725  }
726  //ElementAccessor<3> ele;
727 
728  // modify side fluxes in parallel
729  // for every local edge take time term on digonal and add it to the corresponding flux
730  /*
731  for (unsigned int i_loc = 0; i_loc < el_ds->lsize(); i_loc++) {
732  ele = mesh_->element_accessor(el_4_loc[i_loc]);
733  for (unsigned int i=0; i<ele->n_sides(); i++) {
734  side_rows[i] = side_row_4_id[ mh_dh.side_dof( ele_ac.side(i) ) ];
735  values[i] = -1.0 * ele_ac.measure() *
736  data.cross_section.value(ele_ac.centre(), ele_ac.element_accessor()) *
737  data.water_source_density.value(ele_ac.centre(), ele_ac.element_accessor()) /
738  ele_ac.n_sides();
739  }
740  VecSetValues(schur0->get_solution(), ele_ac.n_sides(), side_rows, values, ADD_VALUES);
741  }
742  VecAssemblyBegin(schur0->get_solution());
743  VecAssemblyEnd(schur0->get_solution());
744  */
745 }
746 
747 
749  START_TIMER("Darcy output data");
750 
751  //print_matlab_matrix("matrix_" + std::to_string(time_->step().index()));
752 
753  //time_->view("DARCY"); //time governor information output
754  this->output_object->output();
755 
756 
757  START_TIMER("Darcy balance output");
758  balance_->calculate_cumulative(data_->water_balance_idx, schur0->get_solution());
759  balance_->calculate_instant(data_->water_balance_idx, schur0->get_solution());
760  balance_->output();
761 }
762 
763 
765 {
766  return schur0->get_solution_precision();
767 }
768 
769 
770 // ===========================================================================================
771 //
772 // MATRIX ASSEMBLY - we use abstract assembly routine, where LS Mat/Vec SetValues
773 // are in fact pointers to allocating or filling functions - this is governed by Linsystem roitunes
774 //
775 // =======================================================================================
777 {
778  START_TIMER("DarcyFlowMHy::assembly_mh_matrix");
779 
780  // set auxiliary flag for switchting Dirichlet like BC
781  data_->force_no_neumann_bc = use_steady_assembly_ && (nonlinear_iteration_ == 0);
782  data_->n_schur_compls = n_schur_compls;
783 
784 
785  balance_->start_flux_assembly(data_->water_balance_idx);
786 
787  // TODO: try to move this into balance, or have it in the generic assembler class, that should perform the cell loop
788  // including various pre- and post-actions
789  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
790  unsigned int dim = dh_cell.dim();
791  assembler[dim-1]->assemble(dh_cell);
792  }
793 
794 
795  balance_->finish_flux_assembly(data_->water_balance_idx);
796 
797 }
798 
799 
801 {
802  START_TIMER("DarcyFlowMH::allocate_mh_matrix");
803 
804  // set auxiliary flag for switchting Dirichlet like BC
805  data_->n_schur_compls = n_schur_compls;
806  LinSys *ls = schur0;
807 
808  // to make space for second schur complement, max. 10 neighbour edges of one el.
809  double zeros[100000];
810  for(int i=0; i<100000; i++) zeros[i] = 0.0;
811 
812  std::vector<LongIdx> tmp_rows;
813  tmp_rows.reserve(200);
814 
815  std::vector<LongIdx> dofs, dofs_ngh;
816  dofs.reserve(data_->dh_->max_elem_dofs());
817  dofs_ngh.reserve(data_->dh_->max_elem_dofs());
818 
819  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
820  ElementAccessor<3> ele = dh_cell.elm();
821 
822  const uint ndofs = dh_cell.n_dofs();
823  dofs.resize(ndofs);
824  dh_cell.get_dof_indices(dofs);
825 
826  // whole local MH matrix
827  ls->mat_set_values(ndofs, dofs.data(), ndofs, dofs.data(), zeros);
828 
829  tmp_rows.clear();
830 
831  // compatible neighborings rows
832  unsigned int n_neighs = ele->n_neighs_vb();
833  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
834  // every compatible connection adds a 2x2 matrix involving
835  // current element pressure and a connected edge pressure
836 
837  // read neighbor dofs (dh dofhandler)
838  // neighbor cell owning neighb_side
839  DHCellAccessor dh_neighb_cell = neighb_side.cell();
840 
841  const uint ndofs_ngh = dh_neighb_cell.n_dofs();
842  dofs_ngh.resize(ndofs_ngh);
843  dh_neighb_cell.get_dof_indices(dofs_ngh);
844 
845  // local index of pedge dof on neighboring cell
846  // (dim+1) is number of edges of higher dim element
847  // TODO: replace with DHCell getter when available for FESystem component
848  const unsigned int t = dh_neighb_cell.n_dofs() - (dh_neighb_cell.dim()+1) + neighb_side.side().side_idx();
849  tmp_rows.push_back(dofs_ngh[t]);
850  }
851 
852  const uint nsides = ele->n_sides();
853  LongIdx * edge_rows = dofs.data() + nsides; // pointer to start of ele
854  // allocate always also for schur 2
855  ls->mat_set_values(nsides+1, edge_rows, n_neighs, tmp_rows.data(), zeros); // (edges, ele) x (neigh edges)
856  ls->mat_set_values(n_neighs, tmp_rows.data(), nsides+1, edge_rows, zeros); // (neigh edges) x (edges, ele)
857  ls->mat_set_values(n_neighs, tmp_rows.data(), n_neighs, tmp_rows.data(), zeros); // (neigh edges) x (neigh edges)
858 
859  tmp_rows.clear();
860 
861  if (data_->mortar_method_ != NoMortar) {
862  auto &isec_list = mesh_->mixed_intersections().element_intersections_[ele.idx()];
863  for(auto &isec : isec_list ) {
864  IntersectionLocalBase *local = isec.second;
865  DHCellAccessor dh_cell_slave = data_->dh_->cell_accessor_from_element(local->bulk_ele_idx());
866 
867  const uint ndofs_slave = dh_cell_slave.n_dofs();
868  dofs_ngh.resize(ndofs_slave);
869  dh_cell_slave.get_dof_indices(dofs_ngh);
870 
871  //DebugOut().fmt("Alloc: {} {}", ele.idx(), local->bulk_ele_idx());
872  for(unsigned int i_side=0; i_side < dh_cell_slave.elm()->n_sides(); i_side++) {
873  // TODO: replace with DHCell getter when available for FESystem component
874  tmp_rows.push_back( dofs_ngh[(ndofs_slave+1)/2+i_side] );
875  //DebugOut() << "aedge" << print_var(tmp_rows[tmp_rows.size()-1]);
876  }
877  }
878  }
879  /*
880  for(unsigned int i_side=0; i_side < ele->n_sides(); i_side++) {
881  DebugOut() << "aedge:" << print_var(edge_rows[i_side]);
882  }*/
883 
884  edge_rows = dofs.data() + nsides +1; // pointer to start of edges
885  ls->mat_set_values(nsides, edge_rows, tmp_rows.size(), tmp_rows.data(), zeros); // master edges x neigh edges
886  ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), nsides, edge_rows, zeros); // neigh edges x master edges
887  ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), tmp_rows.size(), tmp_rows.data(), zeros); // neigh edges x neigh edges
888 
889  }
890 /*
891  // alloc edge diagonal entries
892  if(rank == 0)
893  for( vector<Edge>::iterator edg = mesh_->edges.begin(); edg != mesh_->edges.end(); ++edg) {
894  int edg_idx = mh_dh.row_4_edge[edg->side(0)->edge_idx()];
895 
896 // for( vector<Edge>::iterator edg2 = mesh_->edges.begin(); edg2 != mesh_->edges.end(); ++edg2){
897 // int edg_idx2 = mh_dh.row_4_edge[edg2->side(0)->edge_idx()];
898 // if(edg_idx == edg_idx2){
899 // DBGCOUT(<< "P[ " << rank << " ] " << "edg alloc: " << edg_idx << " " << edg_idx2 << "\n");
900  ls->mat_set_value(edg_idx, edg_idx, 0.0);
901 // }
902 // }
903  }
904  */
905  /*
906  if (mortar_method_ == MortarP0) {
907  P0_CouplingAssembler(*this).assembly(*ls);
908  } else if (mortar_method_ == MortarP1) {
909  P1_CouplingAssembler(*this).assembly(*ls);
910  }*/
911 }
912 
914 {
915  START_TIMER("assembly source term");
916  balance_->start_source_assembly(data_->water_balance_idx);
917 
918  std::vector<LongIdx> global_dofs(data_->dh_->max_elem_dofs());
919 
920  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
921  ElementAccessor<3> ele = dh_cell.elm();
922 
923  const uint ndofs = dh_cell.n_dofs();
924  global_dofs.resize(ndofs);
925  dh_cell.get_dof_indices(global_dofs);
926 
927  double cs = data_->cross_section.value(ele.centre(), ele);
928 
929  // set sources
930  double source = ele.measure() * cs *
931  (data_->water_source_density.value(ele.centre(), ele)
932  +data_->extra_source.value(ele.centre(), ele));
933  // TODO: replace with DHCell getter when available for FESystem component
934  schur0->rhs_set_value(global_dofs[ndofs/2], -1.0 * source );
935 
936  balance_->add_source_values(data_->water_balance_idx, ele.region().bulk_idx(),
937  {dh_cell.get_loc_dof_indices()[ndofs/2]}, {0}, {source});
938  }
939 
940  balance_->finish_source_assembly(data_->water_balance_idx);
941 }
942 
943 
944 
945 
946 /*******************************************************************************
947  * COMPOSE WATER MH MATRIX WITHOUT SCHUR COMPLEMENT
948  ******************************************************************************/
949 
951 
952  START_TIMER("preallocation");
953 
954  if (schur0 == NULL) { // create Linear System for MH matrix
955 
956  if (in_rec.type() == LinSys_BDDC::get_input_type()) {
957 #ifdef FLOW123D_HAVE_BDDCML
958  WarningOut() << "For BDDC no Schur complements are used.";
959  n_schur_compls = 0;
960  LinSys_BDDC *ls = new LinSys_BDDC(&(*data_->dh_->distr()),
961  true); // swap signs of matrix and rhs to make the matrix SPD
962  ls->set_from_input(in_rec);
963  ls->set_solution( data_->full_solution.petsc_vec() );
964  // possible initialization particular to BDDC
965  START_TIMER("BDDC set mesh data");
967  schur0=ls;
968  END_TIMER("BDDC set mesh data");
969 #else
970  Exception
971  xprintf(Err, "Flow123d was not build with BDDCML support.\n");
972 #endif // FLOW123D_HAVE_BDDCML
973  }
974  else if (in_rec.type() == LinSys_PETSC::get_input_type()) {
975  // use PETSC for serial case even when user wants BDDC
976  if (n_schur_compls > 2) {
977  WarningOut() << "Invalid number of Schur Complements. Using 2.";
978  n_schur_compls = 2;
979  }
980 
981  LinSys_PETSC *schur1, *schur2;
982 
983  if (n_schur_compls == 0) {
984  LinSys_PETSC *ls = new LinSys_PETSC( &(*data_->dh_->distr()) );
985 
986  // temporary solution; we have to set precision also for sequantial case of BDDC
987  // final solution should be probably call of direct solver for oneproc case
988  if (in_rec.type() != LinSys_BDDC::get_input_type()) ls->set_from_input(in_rec);
989  else {
990  ls->LinSys::set_from_input(in_rec); // get only common options
991  }
992 
993  ls->set_solution( data_->full_solution.petsc_vec() );
994  schur0=ls;
995  } else {
996  IS is;
997  auto side_dofs_vec = get_component_indices_vec(0);
998 
999  ISCreateGeneral(PETSC_COMM_SELF, side_dofs_vec.size(), &(side_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1000  //ISView(is, PETSC_VIEWER_STDOUT_SELF);
1001  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
1002 
1003  SchurComplement *ls = new SchurComplement(&(*data_->dh_->distr()), is);
1004 
1005  // make schur1
1007  if (n_schur_compls==1) {
1008  schur1 = new LinSys_PETSC(ds);
1009  schur1->set_positive_definite();
1010  } else {
1011  IS is;
1012  auto elem_dofs_vec = get_component_indices_vec(1);
1013 
1014  const PetscInt *b_indices;
1015  ISGetIndices(ls->IsB, &b_indices);
1016  uint b_size = ls->loc_size_B;
1017  for(uint i_b=0, i_bb=0; i_b < b_size && i_bb < elem_dofs_vec.size(); i_b++) {
1018  if (b_indices[i_b] == elem_dofs_vec[i_bb])
1019  elem_dofs_vec[i_bb++] = i_b + ds->begin();
1020  }
1021  ISRestoreIndices(ls->IsB, &b_indices);
1022 
1023 
1024  ISCreateGeneral(PETSC_COMM_SELF, elem_dofs_vec.size(), &(elem_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1025  //ISView(is, PETSC_VIEWER_STDOUT_SELF);
1026  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
1027  SchurComplement *ls1 = new SchurComplement(ds, is); // is is deallocated by SchurComplement
1028  ls1->set_negative_definite();
1029 
1030  // make schur2
1031  schur2 = new LinSys_PETSC( ls1->make_complement_distribution() );
1032  schur2->set_positive_definite();
1033  ls1->set_complement( schur2 );
1034  schur1 = ls1;
1035  }
1036  ls->set_complement( schur1 );
1037  ls->set_from_input(in_rec);
1038  ls->set_solution( data_->full_solution.petsc_vec() );
1039  schur0=ls;
1040  }
1041 
1042  START_TIMER("PETSC PREALLOCATION");
1043  schur0->set_symmetric();
1045 
1047 
1048  VecZeroEntries(schur0->get_solution());
1049  END_TIMER("PETSC PREALLOCATION");
1050  }
1051  else {
1052  xprintf(Err, "Unknown solver type. Internal error.\n");
1053  }
1054 
1055  END_TIMER("preallocation");
1056  }
1057 
1058 }
1059 
1060 
1061 
1062 
1064  START_TIMER("DarcyFlowMH_Steady::assembly_linear_system");
1065 
1066  data_->is_linear=true;
1067  bool is_steady = zero_time_term();
1068  //DebugOut() << "Assembly linear system\n";
1069  if (data_changed_) {
1070  data_changed_ = false;
1071  //DebugOut() << "Data changed\n";
1072  // currently we have no optimization for cases when just time term data or RHS data are changed
1073  START_TIMER("full assembly");
1074  if (typeid(*schur0) != typeid(LinSys_BDDC)) {
1075  schur0->start_add_assembly(); // finish allocation and create matrix
1076  }
1077 
1080 
1082 
1083 
1084  assembly_mh_matrix( data_->multidim_assembler ); // fill matrix
1085 
1087 // print_matlab_matrix("matrix");
1089  //MatView( *const_cast<Mat*>(schur0->get_matrix()), PETSC_VIEWER_STDOUT_WORLD );
1090  //VecView( *const_cast<Vec*>(schur0->get_rhs()), PETSC_VIEWER_STDOUT_WORLD);
1091 
1092  if (! is_steady) {
1093  START_TIMER("fix time term");
1094  //DebugOut() << "setup time term\n";
1095  // assembly time term and rhs
1096  setup_time_term();
1097  modify_system();
1098  }
1099  else
1100  {
1101  balance_->start_mass_assembly(data_->water_balance_idx);
1102  balance_->finish_mass_assembly(data_->water_balance_idx);
1103  }
1104  END_TIMER("full assembly");
1105  } else {
1106  START_TIMER("modify system");
1107  if (! is_steady) {
1108  modify_system();
1109  } else {
1110  //xprintf(PrgErr, "Planned computation time for steady solver, but data are not changed.\n");
1111  }
1112  END_TIMER("modiffy system");
1113  }
1114 
1115 }
1116 
1117 
1118 void DarcyMH::print_matlab_matrix(std::string matlab_file)
1119 {
1120  std::string output_file;
1121 
1122  if ( typeid(*schur0) == typeid(LinSys_BDDC) ){
1123 // WarningOut() << "Can output matrix only on a single processor.";
1124 // output_file = FilePath(matlab_file + "_bddc.m", FilePath::output_file);
1125 // ofstream os( output_file );
1126 // auto bddc = static_cast<LinSys_BDDC*>(schur0);
1127 // bddc->print_matrix(os);
1128  }
1129  else {//if ( typeid(*schur0) == typeid(LinSys_PETSC) ){
1130  output_file = FilePath(matlab_file + ".m", FilePath::output_file);
1131  PetscViewer viewer;
1132  PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
1133  PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1134  MatView( *const_cast<Mat*>(schur0->get_matrix()), viewer);
1135  VecView( *const_cast<Vec*>(schur0->get_rhs()), viewer);
1136  VecView( *const_cast<Vec*>(schur0->get_rhs()), viewer);
1137  VecView( *const_cast<Vec*>(&(schur0->get_solution())), viewer);
1138  }
1139 // else{
1140 // WarningOut() << "No matrix output available for the current solver.";
1141 // return;
1142 // }
1143 
1144  // compute h_min for different dimensions
1145  double d_max = std::numeric_limits<double>::max();
1146  double h1 = d_max, h2 = d_max, h3 = d_max;
1147  double he2 = d_max, he3 = d_max;
1148  for (auto ele : mesh_->elements_range()) {
1149  switch(ele->dim()){
1150  case 1: h1 = std::min(h1,ele.measure()); break;
1151  case 2: h2 = std::min(h2,ele.measure()); break;
1152  case 3: h3 = std::min(h3,ele.measure()); break;
1153  }
1154 
1155  for (unsigned int j=0; j<ele->n_sides(); j++) {
1156  switch(ele->dim()){
1157  case 2: he2 = std::min(he2, ele.side(j)->measure()); break;
1158  case 3: he3 = std::min(he3, ele.side(j)->measure()); break;
1159  }
1160  }
1161  }
1162  if(h1 == d_max) h1 = 0;
1163  if(h2 == d_max) h2 = 0;
1164  if(h3 == d_max) h3 = 0;
1165  if(he2 == d_max) he2 = 0;
1166  if(he3 == d_max) he3 = 0;
1167 
1168  FILE * file;
1169  file = fopen(output_file.c_str(),"a");
1170  fprintf(file, "nA = %d;\n", data_->dh_cr_disc_->distr()->size());
1171  fprintf(file, "nB = %d;\n", data_->dh_->mesh()->get_el_ds()->size());
1172  fprintf(file, "nBF = %d;\n", data_->dh_cr_->distr()->size());
1173  fprintf(file, "h1 = %e;\nh2 = %e;\nh3 = %e;\n", h1, h2, h3);
1174  fprintf(file, "he2 = %e;\nhe3 = %e;\n", he2, he3);
1175  fclose(file);
1176 }
1177 
1178 
1179 //template <int dim>
1180 //std::vector<arma::vec3> dof_points(DHCellAccessor cell, const Mapping<dim, 3> &mapping) {
1181 //
1182 //
1183 // vector<arma::vec::fixed<dim+1>> bary_dof_points = cell->fe()->dof_points();
1184 //
1185 // std::vector<arma::vec3> points(20);
1186 // points.resize(0);
1187 //
1188 //}
1189 //
1190 
1192  START_TIMER("DarcyFlowMH_Steady::set_mesh_data_for_bddc");
1193  // prepare mesh for BDDCML
1194  // initialize arrays
1195  // auxiliary map for creating coordinates of local dofs and global-to-local numbering
1196  std::map<int, arma::vec3> localDofMap;
1197  // connectivity for the subdomain, i.e. global dof numbers on element, stored element-by-element
1198  // Indices of Nodes on Elements
1199  std::vector<int> inet;
1200  // number of degrees of freedom on elements - determines elementwise chunks of INET array
1201  // Number of Nodes on Elements
1202  std::vector<int> nnet;
1203  // Indices of Subdomain Elements in Global Numbering - for local elements, their global indices
1204  std::vector<int> isegn;
1205 
1206  // This array is currently not used in BDDCML, it was used as an interface scaling alternative to scaling
1207  // by diagonal. It corresponds to the rho-scaling.
1208  std::vector<double> element_permeability;
1209 
1210  // maximal and minimal dimension of elements
1211  uint elDimMax = 1;
1212  uint elDimMin = 3;
1213  std::vector<LongIdx> cell_dofs_global(10);
1214 
1215 
1216 
1217  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1218  // for each element, create local numbering of dofs as fluxes (sides), pressure (element centre), Lagrange multipliers (edges), compatible connections
1219 
1220  dh_cell.get_dof_indices(cell_dofs_global);
1221 
1222  inet.insert(inet.end(), cell_dofs_global.begin(), cell_dofs_global.end());
1223  uint n_inet = cell_dofs_global.size();
1224 
1225 
1226  uint dim = dh_cell.elm().dim();
1227  elDimMax = std::max( elDimMax, dim );
1228  elDimMin = std::min( elDimMin, dim );
1229 
1230  // TODO: this is consistent with previous implementation, but may be wrong as it use global element numbering
1231  // used in sequential mesh, do global numbering of distributed elements.
1232  isegn.push_back( dh_cell.elm_idx());
1233 
1234  // TODO: use FiniteElement::dof_points
1235  for (unsigned int si=0; si<dh_cell.elm()->n_sides(); si++) {
1236  arma::vec3 coord = dh_cell.elm().side(si)->centre();
1237  // TODO: replace with DHCell getter when available for FESystem component
1238  // flux dof points
1239  localDofMap.insert( std::make_pair( cell_dofs_global[si], coord ) );
1240  // pressure trace dof points
1241  localDofMap.insert( std::make_pair( cell_dofs_global[si+dim+2], coord ) );
1242  }
1243  // pressure dof points
1244  arma::vec3 elm_centre = dh_cell.elm().centre();
1245  localDofMap.insert( std::make_pair( cell_dofs_global[dim+1], elm_centre ) );
1246 
1247  // insert dofs related to compatible connections
1248  //const Element *ele = dh_cell.elm().element();
1249  for(DHCellSide side : dh_cell.neighb_sides()) {
1250  uint neigh_dim = side.cell().elm().dim();
1251  side.cell().get_dof_indices(cell_dofs_global);
1252  int edge_row = cell_dofs_global[neigh_dim+2+side.side_idx()];
1253  localDofMap.insert( std::make_pair( edge_row, side.centre() ) );
1254  inet.push_back( edge_row );
1255  n_inet++;
1256  }
1257  nnet.push_back(n_inet);
1258 
1259 
1260  // version for rho scaling
1261  // trace computation
1262  double conduct = data_->conductivity.value( elm_centre , dh_cell.elm() );
1263  auto aniso = data_->anisotropy.value( elm_centre , dh_cell.elm() );
1264 
1265  // compute mean on the diagonal
1266  double coef = 0.;
1267  for ( int i = 0; i < 3; i++) {
1268  coef = coef + aniso.at(i,i);
1269  }
1270  // Maybe divide by cs
1271  coef = conduct*coef / 3;
1272 
1273  OLD_ASSERT( coef > 0.,
1274  "Zero coefficient of hydrodynamic resistance %f . \n ", coef );
1275  element_permeability.push_back( 1. / coef );
1276  }
1277 // uint i_inet = 0;
1278 // for(int n_dofs : nnet) {
1279 // DebugOut() << "nnet: " << n_dofs;
1280 // for(int j=0; j < n_dofs; j++, i_inet++) {
1281 // DebugOut() << "inet: " << inet[i_inet];
1282 // }
1283 // }
1284 
1285  auto distr = data_->dh_->distr();
1286 // for(auto pair : localDofMap) {
1287 // DebugOut().every_proc() << "r: " << distr->myp() << " gi: " << pair.first << "xyz: " << pair.second[0];
1288 //
1289 // }
1290 
1291 
1292  //convert set of dofs to vectors
1293  // number of nodes (= dofs) on the subdomain
1294  int numNodeSub = localDofMap.size();
1295  //ASSERT_EQ( (unsigned int)numNodeSub, data_->dh_->lsize() );
1296  // Indices of Subdomain Nodes in Global Numbering - for local nodes, their global indices
1297  std::vector<int> isngn( numNodeSub );
1298  // pseudo-coordinates of local nodes (i.e. dofs)
1299  // they need not be exact, they are used just for some geometrical considerations in BDDCML,
1300  // such as selection of corners maximizing area of a triangle, bounding boxes fro subdomains to
1301  // find candidate neighbours etc.
1302  std::vector<double> xyz( numNodeSub * 3 ) ;
1303  int ind = 0;
1304  std::map<int,arma::vec3>::iterator itB = localDofMap.begin();
1305  for ( ; itB != localDofMap.end(); ++itB ) {
1306  isngn[ind] = itB -> first;
1307 
1308  arma::vec3 coord = itB -> second;
1309  for ( int j = 0; j < 3; j++ ) {
1310  xyz[ j*numNodeSub + ind ] = coord[j];
1311  }
1312 
1313  ind++;
1314  }
1315  localDofMap.clear();
1316 
1317  // Number of Nodal Degrees of Freedom
1318  // nndf is trivially one - dofs coincide with nodes
1319  std::vector<int> nndf( numNodeSub, 1 );
1320 
1321  // prepare auxiliary map for renumbering nodes
1322  typedef std::map<int,int> Global2LocalMap_; //! type for storage of global to local map
1323  Global2LocalMap_ global2LocalNodeMap;
1324  for ( unsigned ind = 0; ind < isngn.size(); ++ind ) {
1325  global2LocalNodeMap.insert( std::make_pair( static_cast<unsigned>( isngn[ind] ), ind ) );
1326  }
1327 
1328  // renumber nodes in the inet array to locals
1329  int indInet = 0;
1330  for ( unsigned int iEle = 0; iEle < isegn.size(); iEle++ ) {
1331  int nne = nnet[ iEle ];
1332  for ( int ien = 0; ien < nne; ien++ ) {
1333 
1334  int indGlob = inet[indInet];
1335  // map it to local node
1336  Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1337  OLD_ASSERT( pos != global2LocalNodeMap.end(),
1338  "Cannot remap node index %d to local indices. \n ", indGlob );
1339  int indLoc = static_cast<int> ( pos -> second );
1340 
1341  // store the node
1342  inet[ indInet++ ] = indLoc;
1343  }
1344  }
1345 
1346  int numNodes = size;
1347  int numDofsInt = size;
1348  int spaceDim = 3; // TODO: what is the proper value here?
1349  int meshDim = elDimMax;
1350 
1351  /**
1352  * We need:
1353  * - local to global element map (possibly mesh->el_4_loc
1354  * - inet, nnet - local dof numbers per element, local numbering of only those dofs that are on owned elements
1355  * 1. collect DH local dof indices on elements, manage map from DH local indices to BDDC local dof indices
1356  * 2. map collected DH indices to BDDC indices using the map
1357  * - local BDDC dofs to global dofs, use DH to BDDC map with DH local to global map
1358  * - XYZ - permuted, collect in main loop into array of size of all DH local dofs, compress and rearrange latter
1359  * - element_permeability - in main loop
1360  */
1361  bddc_ls -> load_mesh( LinSys_BDDC::BDDCMatrixType::SPD_VIA_SYMMETRICGENERAL, spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1362 }
1363 
1364 
1365 
1366 
1367 //=============================================================================
1368 // DESTROY WATER MH SYSTEM STRUCTURE
1369 //=============================================================================
1371  if (previous_solution != nullptr) chkerr(VecDestroy(&previous_solution));
1372  if (steady_diagonal != nullptr) chkerr(VecDestroy(&steady_diagonal));
1373  if (new_diagonal != nullptr) chkerr(VecDestroy(&new_diagonal));
1374  if (steady_rhs != nullptr) chkerr(VecDestroy(&steady_rhs));
1375 
1376 
1377  if (schur0 != NULL) {
1378  delete schur0;
1379  }
1380 
1381  if (output_object) delete output_object;
1382 
1383  if(time_ != nullptr)
1384  delete time_;
1385 
1386 }
1387 
1388 
1389 /*
1390 void mat_count_off_proc_values(Mat m, Vec v) {
1391  int n, first, last;
1392  const PetscInt *cols;
1393  Distribution distr(v);
1394 
1395  int n_off = 0;
1396  int n_on = 0;
1397  int n_off_rows = 0;
1398  MatGetOwnershipRange(m, &first, &last);
1399  for (int row = first; row < last; row++) {
1400  MatGetRow(m, row, &n, &cols, PETSC_NULL);
1401  bool exists_off = false;
1402  for (int i = 0; i < n; i++)
1403  if (distr.get_proc(cols[i]) != distr.myp())
1404  n_off++, exists_off = true;
1405  else
1406  n_on++;
1407  if (exists_off)
1408  n_off_rows++;
1409  MatRestoreRow(m, row, &n, &cols, PETSC_NULL);
1410  }
1411 }
1412 */
1413 
1414 
1415 
1416 
1417 
1418 
1419 
1420 
1421 
1422 
1423 
1424 
1426 {
1427  double *local_sol = schur0->get_solution_array();
1428 
1429  // cycle over local element rows
1430 
1431  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1432  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1433  ElementAccessor<3> ele = dh_cell.elm();
1434  // set initial condition
1435  // TODO: replace with DHCell getter when available for FESystem component
1436  const IntIdx p_ele_dof = dh_cell.get_loc_dof_indices()[dh_cell.n_dofs()/2];
1437  local_sol[p_ele_dof] = data_->init_pressure.value(ele.centre(),ele);
1438  }
1439 }
1440 
1442  // save diagonal of steady matrix
1443  MatGetDiagonal(*( schur0->get_matrix() ), steady_diagonal);
1444  // save RHS
1445  VecCopy(*( schur0->get_rhs()), steady_rhs);
1446 
1447 
1448  PetscScalar *local_diagonal;
1449  VecGetArray(new_diagonal,& local_diagonal);
1450 
1451  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1452 
1453  balance_->start_mass_assembly(data_->water_balance_idx);
1454 
1455  std::vector<LongIdx> dofs;
1456  dofs.reserve(data_->dh_->max_elem_dofs());
1457 
1458  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1459  ElementAccessor<3> ele = dh_cell.elm();
1460  dofs.resize(dh_cell.n_dofs());
1461  dh_cell.get_dof_indices(dofs);
1462 
1463  // TODO: replace with DHCell getter when available for FESystem component
1464  const uint p_ele_dof = dh_cell.n_dofs() / 2;
1465  // set new diagonal
1466  double diagonal_coeff = data_->cross_section.value(ele.centre(), ele)
1467  * ( data_->storativity.value(ele.centre(), ele)
1468  +data_->extra_storativity.value(ele.centre(), ele))
1469  * ele.measure();
1470  local_diagonal[dh_cell.get_loc_dof_indices()[p_ele_dof]]= - diagonal_coeff / time_->dt();
1471 
1472  balance_->add_mass_values(data_->water_balance_idx, dh_cell,
1473  {dh_cell.get_loc_dof_indices()[p_ele_dof]},
1474  {diagonal_coeff}, 0.0);
1475  }
1476  VecRestoreArray(new_diagonal,& local_diagonal);
1477  MatDiagonalSet(*( schur0->get_matrix() ), new_diagonal, ADD_VALUES);
1478 
1480 
1481  balance_->finish_mass_assembly(data_->water_balance_idx);
1482 }
1483 
1485  START_TIMER("modify system");
1486  if (time_->is_changed_dt() && time_->step().index()>0) {
1487  double scale_factor=time_->step(-2).length()/time_->step().length();
1488  if (scale_factor != 1.0) {
1489  // if time step has changed and setup_time_term not called
1490  MatDiagonalSet(*( schur0->get_matrix() ),steady_diagonal, INSERT_VALUES);
1491 
1492  VecScale(new_diagonal, time_->last_dt()/time_->dt());
1493  MatDiagonalSet(*( schur0->get_matrix() ),new_diagonal, ADD_VALUES);
1495  }
1496  }
1497 
1498  // modify RHS - add previous solution
1499  VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, previous_solution);
1500 // VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, schur0->get_solution());
1501  VecAXPY(*( schur0->get_rhs()), 1.0, steady_rhs);
1503 
1504  //VecSwap(previous_solution, schur0->get_solution());
1505 }
1506 
1507 
1508 /// Helper method fills range (min and max) of given component
1509 void dofs_range(unsigned int n_dofs, unsigned int &min, unsigned int &max, unsigned int component) {
1510  if (component==0) {
1511  min = 0;
1512  max = n_dofs/2;
1513  } else if (component==1) {
1514  min = n_dofs/2;
1515  max = (n_dofs+1)/2;
1516  } else {
1517  min = (n_dofs+1)/2;
1518  max = n_dofs;
1519  }
1520 }
1521 
1522 
1524  ASSERT_LT_DBG(component, 3).error("Invalid component!");
1525  unsigned int i, n_dofs, min, max;
1526  std::vector<int> dof_vec;
1527  std::vector<LongIdx> dof_indices(data_->dh_->max_elem_dofs());
1528  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1529  n_dofs = dh_cell.get_dof_indices(dof_indices);
1530  dofs_range(n_dofs, min, max, component);
1531  for (i=min; i<max; ++i) dof_vec.push_back(dof_indices[i]);
1532  }
1533  return dof_vec;
1534 }
1535 
1536 
1537 //-----------------------------------------------------------------------------
1538 // vim: set cindent:
TimeGovernor & time()
Definition: equation.hh:149
Field< 3, FieldValue< 3 >::VectorFixed > field_ele_velocity
Field< 3, FieldValue< 3 >::Scalar > extra_source
Externally added storativity.
FieldSet * eq_data_
Definition: equation.hh:227
std::vector< int > get_component_indices_vec(unsigned int component) const
Get vector of all DOF indices of given component (0..side, 1..element, 2..edge)
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:50
Output class for darcy_flow_mh model.
Field< 3, FieldValue< 3 >::Scalar > init_pressure
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:328
virtual void initialize_specific()
void set_complement(LinSys_PETSC *ls)
Set complement LinSys object.
Definition: schur.cc:287
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
ArmaVec< double, N > vec
Definition: armor.hh:861
SchurComplement SchurComplement
Definition: linsys.hh:92
void assembly_mh_matrix(MultidimAssembly &assembler)
void set_symmetric(bool flag=true)
Definition: linsys.hh:561
unsigned int uint
Classes with algorithms for computation of intersections of meshes.
FieldCommon & input_selection(Input::Type::Selection element_selection)
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:602
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:822
static const Input::Type::Record & get_input_type()
The specification of output stream.
Definition: output_time.cc:37
Common base for intersection object.
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
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:452
void set_from_input(const Input::Record in_rec) override
int IntIdx
Definition: index_types.hh:25
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:267
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
Abstract linear system class.
Definition: balance.hh:40
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_
Field< 3, FieldValue< 3 >::Scalar > extra_storativity
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
Field< 3, FieldValue< 3 >::Scalar > sigma
unsigned int dim() const
Return dimension of element appropriate to cell.
Vec steady_rhs
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
void initialize() override
void set_negative_definite(bool flag=true)
Definition: linsys.hh:588
Definition: mesh.h:78
static const Input::Type::Record & type_field_descriptor()
virtual void start_allocation()
Definition: linsys.hh:333
Cell accessor allow iterate over DOF handler cells.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:148
virtual void finish_assembly()=0
virtual double get_solution_precision()=0
std::vector< std::vector< ILpair > > element_intersections_
static const Input::Type::Instance & get_input_type(FieldSet &eq_data, const std::string &equation_name)
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:347
const TimeStep & step(int index=-1) const
static const std::string field_descriptor_record_description(const string &record_name)
Definition: field_common.cc:73
int loc_size_B
Definition: schur.hh:137
Class for declaration of the integral input data.
Definition: type_base.hh:483
Basic time management functionality for unsteady (and steady) solvers (class Equation).
unsigned int bulk_ele_idx() const
Returns index of bulk element.
unsigned int n_dofs() const
Return number of dofs on given cell.
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
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:38
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.
std::shared_ptr< DiscreteSpace > ds
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
Field< 3, FieldValue< 3 >::Scalar > field_ele_pressure
Externally added water source.
void view(const char *name="") const
unsigned int min_n_it_
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Definition: equation.cc:35
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...
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
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
void solve_time_step(bool output=true)
Solve the problem without moving to next time and without output.
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void set_from_input(const Input::Record in_rec) override
Definition: linsys_BDDC.cc:278
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
virtual void postprocess()
postprocess velocity field (add sources)
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
Field< 3, FieldValue< 3 >::Scalar > water_source_density
static const Input::Type::Instance & get_input_type_specific()
double * get_solution_array()
Definition: linsys.hh:325
BCField< 3, FieldValue< 3 >::Scalar > bc_switch_pressure
Field< 3, FieldValue< 3 >::Scalar > storativity
static constexpr Mask input_copy
Definition: field_flag.hh:44
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
std::shared_ptr< Balance > balance_
#define MPI_MIN
Definition: mpi.h:198
int n_schur_compls
FieldCommon & input_default(const string &input_default)
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
Definitions of basic Lagrangean finite elements with polynomial shape functions.
unsigned int begin(int proc) const
get starting local index
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
const Ret val(const string &key) const
unsigned int n_sides() const
Definition: elements.h:132
#define xprintf(...)
Definition: system.hh:93
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
void set_from_input(const Input::Record in_rec) override
Definition: schur.cc:138
Mesh * mesh_
Definition: equation.hh:218
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.
bool data_changed_
double length() 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:1153
bool use_steady_assembly_
double measure() const
Computes the measure of the element.
unsigned int n_sides() const
Definition: mesh.cc:230
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:458
Region region() const
Definition: accessors.hh:165
LinSys * schur0
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:273
void set_solution(Vec sol_vec)
Definition: linsys.hh:290
Field< 3, FieldValue< 3 >::VectorFixed > flux
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
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.
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
FieldCommon & description(const string &description)
void create_linear_system(Input::AbstractRecord rec)
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:344
Distribution * make_complement_distribution()
get distribution of complement object if complement is defined
Definition: schur.cc:294
const Selection & close() const
Close the Selection, no more values can be added.
Input::Record input_record_
Definition: equation.hh:220
auto disable_where(const Field< spacedim, typename FieldValue< spacedim >::Enum > &control_field, const vector< FieldEnum > &value_list) -> Field &
Definition: field.impl.hh:178
#define MPI_INT
Definition: mpi.h:160
double dt() const
Vec steady_diagonal
Field< 3, FieldValue< 3 >::TensorFixed > anisotropy
DarcyMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:270
FieldCommon & name(const string &name)
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:65
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
Field< 3, FieldValue< 3 >::Scalar > field_edge_pressure
void zero_time_step() override
void dofs_range(unsigned int n_dofs, unsigned int &min, unsigned int &max, unsigned int component)
Helper method fills range (min and max) of given component.
void set_positive_definite(bool flag=true)
Definition: linsys.hh:576
virtual void read_initial_condition()
Record type proxy class.
Definition: type_record.hh:182
FieldCommon & set_limits(double min, double max=std::numeric_limits< double >::max())
Field< 3, FieldValue< 3 >::Scalar > cross_section
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
BCField< 3, FieldValue< 3 >::Enum > bc_type
DarcyFlowMHOutput * output_object
Field< 3, FieldValue< 3 >::Scalar > conductivity
FieldCommon & flags(FieldFlag::Flags::Mask mask)
unsigned int n_edges() const
Definition: mesh.h:134
void prepare_new_time_step()
Class for representation SI units of Fields.
Definition: unit_si.hh:40
virtual const Mat * get_matrix()
Definition: linsys.hh:187
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
double last_dt() const
Vec new_diagonal
virtual double solution_precision() const
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
unsigned int n_neighs_vb() const
Return number of neighbours.
Definition: elements.h:67
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:276
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
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.
Side accessor allows to iterate over sides of DOF handler cell.
MortarMethod mortar_method_
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:180
void rhs_set_value(int row, double val)
Definition: linsys.hh:381
double tolerance_
Input::Type::Record make_field_descriptor_type(const std::string &equation_name) const
Definition: field_set.cc:62
static const Input::Type::Record & get_input_type()
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
TimeGovernor * time_
Definition: equation.hh:219
virtual bool zero_time_term(bool time_global=false)
void output()
Calculate values for output.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:300
Mixed-hybrid model of linear Darcy flow, possibly unsteady.