Flow123d  JS_before_hm-919-g5f1bbbf
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 += anisotropy.name("anisotropy")
210  .description("Anisotropy of the conductivity tensor.")
211  .input_default("1.0")
213 
214  *this += cross_section.name("cross_section")
215  .description("Complement dimension parameter (cross section for 1D, thickness for 2D).")
216  .input_default("1.0")
217  .units( UnitSI().m(3).md() );
218 
219  *this += conductivity.name("conductivity")
220  .description("Isotropic conductivity scalar.")
221  .input_default("1.0")
222  .units( UnitSI().m().s(-1) )
223  .set_limits(0.0);
224 
225  *this += sigma.name("sigma")
226  .description("Transition coefficient between dimensions.")
227  .input_default("1.0")
229 
230  *this += water_source_density.name("water_source_density")
231  .description("Water source density.")
232  .input_default("0.0")
233  .units( UnitSI().s(-1) );
234 
235  *this += bc_type.name("bc_type")
236  .description("Boundary condition type.")
238  .input_default("\"none\"")
240 
241  *this += bc_pressure
243  .name("bc_pressure")
244  .description("Prescribed pressure value on the boundary. Used for all values of ``bc_type`` except ``none`` and ``seepage``. "
245  "See documentation of ``bc_type`` for exact meaning of ``bc_pressure`` in individual boundary condition types.")
246  .input_default("0.0")
247  .units( UnitSI().m() );
248 
249  *this += bc_flux
251  .name("bc_flux")
252  .description("Incoming water boundary flux. Used for bc_types : ``total_flux``, ``seepage``, ``river``.")
253  .input_default("0.0")
254  .units( UnitSI().m().s(-1) );
255 
256  *this += bc_robin_sigma
258  .name("bc_robin_sigma")
259  .description("Conductivity coefficient in the ``total_flux`` or the ``river`` boundary condition type.")
260  .input_default("0.0")
261  .units( UnitSI().s(-1) );
262 
263  *this += bc_switch_pressure
265  .name("bc_switch_pressure")
266  .description("Critical switch pressure for ``seepage`` and ``river`` boundary conditions.")
267  .input_default("0.0")
268  .units( UnitSI().m() );
269 
270 
271  //these are for unsteady
272  *this += init_pressure.name("init_pressure")
273  .description("Initial condition for pressure in time dependent problems.")
274  .input_default("0.0")
275  .units( UnitSI().m() );
276 
277  *this += storativity.name("storativity")
278  .description("Storativity (in time dependent problems).")
279  .input_default("0.0")
280  .units( UnitSI().m(-1) );
281 
282  *this += extra_storativity.name("extra_storativity")
283  .description("Storativity added from upstream equation.")
284  .units( UnitSI().m(-1) )
285  .input_default("0.0")
286  .flags( input_copy );
287 
288  *this += extra_source.name("extra_water_source_density")
289  .description("Water source density added from upstream equation.")
290  .input_default("0.0")
291  .units( UnitSI().s(-1) )
292  .flags( input_copy );
293 
294  //time_term_fields = this->subset({"storativity"});
295  //main_matrix_fields = this->subset({"anisotropy", "conductivity", "cross_section", "sigma", "bc_type", "bc_robin_sigma"});
296  //rhs_fields = this->subset({"water_source_density", "bc_pressure", "bc_flux"});
297 }
298 
299 
300 
301 //=============================================================================
302 // CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
303 // - do it in parallel:
304 // - initial distribution of elements, edges
305 //
306 /*! @brief CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
307  *
308  * Parameters {Solver,NSchurs} number of performed Schur
309  * complements (0,1,2) for water flow MH-system
310  *
311  */
312 //=============================================================================
313 DarcyMH::DarcyMH(Mesh &mesh_in, const Input::Record in_rec, TimeGovernor *tm)
314 : DarcyFlowInterface(mesh_in, in_rec),
315  output_object(nullptr),
316  data_changed_(false),
317  schur0(nullptr),
318  steady_diagonal(nullptr),
319  steady_rhs(nullptr),
320  new_diagonal(nullptr),
321  previous_solution(nullptr)
322 {
323 
324  START_TIMER("Darcy constructor");
325  {
326  auto time_rec = in_rec.val<Input::Record>("time");
327  if (tm == nullptr)
328  {
329  time_ = new TimeGovernor(time_rec);
330  }
331  else
332  {
333  TimeGovernor tm_from_rec(time_rec);
334  if (!tm_from_rec.is_default()) // is_default() == false when time record is present in input file
335  {
336  MessageOut() << "Duplicate key 'time', time in flow equation is already initialized from parent class!";
337  ASSERT(false);
338  }
339  time_ = tm;
340  }
341  }
342 
343  data_ = make_shared<EqData>();
345 
346  data_->is_linear=true;
347 
348  size = mesh_->n_elements() + mesh_->n_sides() + mesh_->n_edges();
349  n_schur_compls = in_rec.val<int>("n_schurs");
350  data_->mortar_method_= in_rec.val<MortarMethod>("mortar_method");
351  if (data_->mortar_method_ != NoMortar) {
353  }
354 
355 
356 
357  //side_ds->view( std::cout );
358  //el_ds->view( std::cout );
359  //edge_ds->view( std::cout );
360  //rows_ds->view( std::cout );
361 
362 }
363 
364 
365 
367 //connecting data fields with mesh
368 {
369 
370  START_TIMER("data init");
371  data_->mesh = mesh_;
372  data_->set_mesh(*mesh_);
373 
374  auto gravity_array = input_record_.val<Input::Array>("gravity");
375  std::vector<double> gvec;
376  gravity_array.copy_to(gvec);
377  gvec.push_back(0.0); // zero pressure shift
378  data_->gravity_ = arma::vec(gvec);
379  data_->gravity_vec_ = data_->gravity_.subvec(0,2);
380 
381  data_->bc_pressure.add_factory(
382  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
383  (data_->gravity_, "bc_piezo_head") );
384  data_->bc_switch_pressure.add_factory(
385  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
386  (data_->gravity_, "bc_switch_piezo_head") );
387  data_->init_pressure.add_factory(
388  std::make_shared<FieldAddPotential<3, FieldValue<3>::Scalar>::FieldFactory>
389  (data_->gravity_, "init_piezo_head") );
390 
391 
392  data_->set_input_list( this->input_record_.val<Input::Array>("input_fields"), *time_ );
393  // Check that the time step was set for the transient simulation.
394  if (! zero_time_term(true) && time_->is_default() ) {
395  //THROW(ExcAssertMsg());
396  //THROW(ExcMissingTimeGovernor() << input_record_.ei_address());
397  MessageOut() << "Missing the key 'time', obligatory for the transient problems." << endl;
398  ASSERT(false);
399  }
400 
401  data_->mark_input_times(*time_);
402 }
403 
404 
405 
407 
408  { // init DOF handler for pressure fields
409 // std::shared_ptr< FiniteElement<0> > fe0_rt = std::make_shared<FE_RT0_disc<0>>();
410  std::shared_ptr< FiniteElement<1> > fe1_rt = std::make_shared<FE_RT0_disc<1>>();
411  std::shared_ptr< FiniteElement<2> > fe2_rt = std::make_shared<FE_RT0_disc<2>>();
412  std::shared_ptr< FiniteElement<3> > fe3_rt = std::make_shared<FE_RT0_disc<3>>();
413  std::shared_ptr< FiniteElement<0> > fe0_disc = std::make_shared<FE_P_disc<0>>(0);
414  std::shared_ptr< FiniteElement<1> > fe1_disc = std::make_shared<FE_P_disc<1>>(0);
415  std::shared_ptr< FiniteElement<2> > fe2_disc = std::make_shared<FE_P_disc<2>>(0);
416  std::shared_ptr< FiniteElement<3> > fe3_disc = std::make_shared<FE_P_disc<3>>(0);
417  std::shared_ptr< FiniteElement<0> > fe0_cr = std::make_shared<FE_CR<0>>();
418  std::shared_ptr< FiniteElement<1> > fe1_cr = std::make_shared<FE_CR<1>>();
419  std::shared_ptr< FiniteElement<2> > fe2_cr = std::make_shared<FE_CR<2>>();
420  std::shared_ptr< FiniteElement<3> > fe3_cr = std::make_shared<FE_CR<3>>();
421 // static FiniteElement<0> fe0_sys = FE_P_disc<0>(0); //TODO fix and use solution with FESystem<0>( {fe0_rt, fe0_disc, fe0_cr} )
422  FESystem<0> fe0_sys( {fe0_disc, fe0_disc, fe0_cr} );
423  FESystem<1> fe1_sys( {fe1_rt, fe1_disc, fe1_cr} );
424  FESystem<2> fe2_sys( {fe2_rt, fe2_disc, fe2_cr} );
425  FESystem<3> fe3_sys( {fe3_rt, fe3_disc, fe3_cr} );
426  MixedPtr<FESystem> fe_sys( std::make_shared<FESystem<0>>(fe0_sys), std::make_shared<FESystem<1>>(fe1_sys),
427  std::make_shared<FESystem<2>>(fe2_sys), std::make_shared<FESystem<3>>(fe3_sys) );
428  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_sys);
429  data_->dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
430  data_->dh_->distribute_dofs(ds);
431  }
432 
433  init_eq_data();
434  this->data_->multidim_assembler = AssemblyBase::create< AssemblyMH >(data_);
436 
437  { // construct pressure, velocity and piezo head fields
438  uint rt_component = 0;
439  ele_flux_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(data_->dh_, rt_component);
440  auto ele_velocity_ptr = std::make_shared< FieldDivide<3, FieldValue<3>::VectorFixed> >(ele_flux_ptr, data_->cross_section);
441  data_->field_ele_velocity.set_field(mesh_->region_db().get_region_set("ALL"), ele_velocity_ptr);
442  data_->full_solution = ele_flux_ptr->get_data_vec();
443 
444  uint p_ele_component = 0;
445  auto ele_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(data_->dh_, p_ele_component, &data_->full_solution);
446  data_->field_ele_pressure.set_field(mesh_->region_db().get_region_set("ALL"), ele_pressure_ptr);
447 
448  uint p_edge_component = 1;
449  auto edge_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(data_->dh_, p_edge_component, &data_->full_solution);
450  data_->field_edge_pressure.set_field(mesh_->region_db().get_region_set("ALL"), edge_pressure_ptr);
451 
452  arma::vec4 gravity = (-1) * data_->gravity_;
453  auto ele_piezo_head_ptr = std::make_shared< FieldAddPotential<3, FieldValue<3>::Scalar> >(gravity, ele_pressure_ptr);
454  data_->field_ele_piezo_head.set_field(mesh_->region_db().get_region_set("ALL"), ele_piezo_head_ptr);
455  }
456 
457  { // init DOF handlers represents edge DOFs
458  uint p_edge_component = 2;
459  data_->dh_cr_ = std::make_shared<SubDOFHandlerMultiDim>(data_->dh_,p_edge_component);
460  }
461 
462  { // init DOF handlers represents side DOFs
463  MixedPtr<FE_CR_disc> fe_cr_disc;
464  std::shared_ptr<DiscreteSpace> ds_cr_disc = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_cr_disc);
465  data_->dh_cr_disc_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
466  data_->dh_cr_disc_->distribute_dofs(ds_cr_disc);
467  }
468 
469  // Initialize bc_switch_dirichlet to size of global boundary.
470  data_->bc_switch_dirichlet.resize(mesh_->n_elements()+mesh_->n_elements(true), 1);
471 
472 
475  .val<Input::Record>("nonlinear_solver")
476  .val<Input::AbstractRecord>("linear_solver");
477 
478  // auxiliary set_time call since allocation assembly evaluates fields as well
481 
482 
483 
484  // allocate time term vectors
485  VecDuplicate(schur0->get_solution(), &previous_solution);
486  VecCreateMPI(PETSC_COMM_WORLD, data_->dh_->distr()->lsize(),PETSC_DETERMINE,&(steady_diagonal));
487  VecDuplicate(steady_diagonal,& new_diagonal);
488  VecZeroEntries(new_diagonal);
489  VecDuplicate(steady_diagonal, &steady_rhs);
490 
491 
492  // initialization of balance object
493  balance_ = std::make_shared<Balance>("water", mesh_);
494  balance_->init_from_input(input_record_.val<Input::Record>("balance"), time());
495  data_->water_balance_idx = balance_->add_quantity("water_volume");
496  balance_->allocate(data_->dh_->distr()->lsize(), 1);
497  balance_->units(UnitSI().m(3));
498 
499 
500  data_->balance = balance_;
501  data_->lin_sys = schur0;
502 
503 
505 }
506 
508 {}
509 
511 {
512 
513  /* TODO:
514  * - Allow solution reconstruction (pressure and velocity) from initial condition on user request.
515  * - Steady solution as an intitial condition may be forced by setting inti_time =-1, and set data for the steady solver in that time.
516  * Solver should be able to switch from and to steady case depending on the zero time term.
517  */
518 
520 
521  // zero_time_term means steady case
522  bool zero_time_term_from_right = zero_time_term();
523 
524 
525  if (zero_time_term_from_right) {
526  // steady case
527  VecZeroEntries(schur0->get_solution());
528  //read_initial_condition(); // Possible solution guess for steady case.
529  use_steady_assembly_ = true;
530  solve_nonlinear(); // with right limit data
531  } else {
532  VecZeroEntries(schur0->get_solution());
533  VecZeroEntries(previous_solution);
535  assembly_linear_system(); // in particular due to balance
536 // print_matlab_matrix("matrix_zero");
537  // TODO: reconstruction of solution in zero time.
538  }
539  //solution_output(T,right_limit); // data for time T in any case
540  output_data();
541 }
542 
543 //=============================================================================
544 // COMPOSE and SOLVE WATER MH System possibly through Schur complements
545 //=============================================================================
547 {
548  START_TIMER("Solving MH system");
549 
550  time_->next_time();
551 
552  time_->view("DARCY"); //time governor information output
553 
554  solve_time_step();
555 }
556 
557 
559 {
561 
562  bool zero_time_term_from_left=zero_time_term();
563 
564  bool jump_time = data_->storativity.is_jump_time();
565  if (! zero_time_term_from_left) {
566  // time term not treated as zero
567  // Unsteady solution up to the T.
568 
569  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
570  use_steady_assembly_ = false;
571 
573  solve_nonlinear(); // with left limit data
574  if (jump_time) {
575  WarningOut() << "Output of solution discontinuous in time not supported yet.\n";
576  //solution_output(T, left_limit); // output use time T- delta*dt
577  //output_data();
578  }
579  }
580 
581  if (time_->is_end()) {
582  // output for unsteady case, end_time should not be the jump time
583  // but rether check that
584  if (! zero_time_term_from_left && ! jump_time && output) output_data();
585  return;
586  }
587 
589  bool zero_time_term_from_right=zero_time_term();
590  if (zero_time_term_from_right) {
591  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
592  use_steady_assembly_ = true;
593  solve_nonlinear(); // with right limit data
594 
595  } else if (! zero_time_term_from_left && jump_time) {
596  WarningOut() << "Discontinuous time term not supported yet.\n";
597  //solution_transfer(); // internally call set_time(T, left) and set_time(T,right) again
598  //solve_nonlinear(); // with right limit data
599  }
600 
601  //solution_output(T,right_limit); // data for time T in any case
602  if (output) output_data();
603 }
604 
605 
606 bool DarcyMH::zero_time_term(bool time_global) {
607  if (time_global) {
608  return (data_->storativity.input_list_size() == 0);
609  } else {
610  return data_->storativity.field_result(mesh_->region_db().get_region_set("BULK")) == result_zeros;
611  }
612 }
613 
614 
616 {
617 
619  double residual_norm = schur0->compute_residual();
621  MessageOut().fmt("[nonlinear solver] norm of initial residual: {}\n", residual_norm);
622 
623  // Reduce is_linear flag.
624  int is_linear_common;
625  MPI_Allreduce(&(data_->is_linear), &is_linear_common,1, MPI_INT ,MPI_MIN,PETSC_COMM_WORLD);
626 
627  Input::Record nl_solver_rec = input_record_.val<Input::Record>("nonlinear_solver");
628  this->tolerance_ = nl_solver_rec.val<double>("tolerance");
629  this->max_n_it_ = nl_solver_rec.val<unsigned int>("max_it");
630  this->min_n_it_ = nl_solver_rec.val<unsigned int>("min_it");
631  if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->max_n_it_;
632 
633  if (! is_linear_common) {
634  // set tolerances of the linear solver unless they are set by user.
635  schur0->set_tolerances(0.1*this->tolerance_, 0.01*this->tolerance_, 100);
636  }
637  vector<double> convergence_history;
638 
639  Vec save_solution;
640  VecDuplicate(schur0->get_solution(), &save_solution);
641  while (nonlinear_iteration_ < this->min_n_it_ ||
642  (residual_norm > this->tolerance_ && nonlinear_iteration_ < this->max_n_it_ )) {
643  OLD_ASSERT_EQUAL( convergence_history.size(), nonlinear_iteration_ );
644  convergence_history.push_back(residual_norm);
645 
646  // stagnation test
647  if (convergence_history.size() >= 5 &&
648  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
649  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
650  // stagnation
651  if (input_record_.val<Input::Record>("nonlinear_solver").val<bool>("converge_on_stagnation")) {
652  WarningOut().fmt("Accept solution on stagnation. Its: {} Residual: {}\n", nonlinear_iteration_, residual_norm);
653  break;
654  } else {
655  THROW(ExcSolverDiverge() << EI_Reason("Stagnation."));
656  }
657  }
658 
659  if (! is_linear_common)
660  VecCopy( schur0->get_solution(), save_solution);
663 
664  // hack to make BDDC work with empty compute_residual
665  if (is_linear_common){
666  // we want to print this info in linear (and steady) case
667  residual_norm = schur0->compute_residual();
668  MessageOut().fmt("[nonlinear solver] lin. it: {}, reason: {}, residual: {}\n",
669  si.n_iterations, si.converged_reason, residual_norm);
670  break;
671  }
672  data_changed_=true; // force reassembly for non-linear case
673 
674  double alpha = 1; // how much of new solution
675  VecAXPBY(schur0->get_solution(), (1-alpha), alpha, save_solution);
676 
677  //LogOut().fmt("Linear solver ended with reason: {} \n", si.converged_reason );
678  //OLD_ASSERT( si.converged_reason >= 0, "Linear solver failed to converge. Convergence reason %d \n", si.converged_reason );
680 
681  residual_norm = schur0->compute_residual();
682  MessageOut().fmt("[nonlinear solver] it: {} lin. it: {}, reason: {}, residual: {}\n",
683  nonlinear_iteration_, si.n_iterations, si.converged_reason, residual_norm);
684  }
685  chkerr(VecDestroy(&save_solution));
686  this -> postprocess();
687 
688  // adapt timestep
689  if (! this->zero_time_term()) {
690  double mult = 1.0;
691  if (nonlinear_iteration_ < 3) mult = 1.6;
692  if (nonlinear_iteration_ > 7) mult = 0.7;
693  time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
694  // int result = time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
695  //DebugOut().fmt("time adaptivity, res: {} it: {} m: {} dt: {} edt: {}\n", result, nonlinear_iteration_, mult, time_->dt(), time_->estimate_dt());
696  }
697 }
698 
700 {
702 }
703 
705 {
706  START_TIMER("postprocess");
707 
708  //fix velocity when mortar method is used
709  if(data_->mortar_method_ != MortarMethod::NoMortar){
710  auto multidim_assembler = AssemblyBase::create< AssemblyMH >(data_);
711  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
712  unsigned int dim = dh_cell.dim();
713  multidim_assembler[dim-1]->fix_velocity(dh_cell);
714  }
715  }
716  //ElementAccessor<3> ele;
717 
718  // modify side fluxes in parallel
719  // for every local edge take time term on digonal and add it to the corresponding flux
720  /*
721  for (unsigned int i_loc = 0; i_loc < el_ds->lsize(); i_loc++) {
722  ele = mesh_->element_accessor(el_4_loc[i_loc]);
723  for (unsigned int i=0; i<ele->n_sides(); i++) {
724  side_rows[i] = side_row_4_id[ mh_dh.side_dof( ele_ac.side(i) ) ];
725  values[i] = -1.0 * ele_ac.measure() *
726  data.cross_section.value(ele_ac.centre(), ele_ac.element_accessor()) *
727  data.water_source_density.value(ele_ac.centre(), ele_ac.element_accessor()) /
728  ele_ac.n_sides();
729  }
730  VecSetValues(schur0->get_solution(), ele_ac.n_sides(), side_rows, values, ADD_VALUES);
731  }
732  VecAssemblyBegin(schur0->get_solution());
733  VecAssemblyEnd(schur0->get_solution());
734  */
735 }
736 
737 
739  START_TIMER("Darcy output data");
740 
741  //print_matlab_matrix("matrix_" + std::to_string(time_->step().index()));
742 
743  //time_->view("DARCY"); //time governor information output
744  this->output_object->output();
745 
746 
747  START_TIMER("Darcy balance output");
748  balance_->calculate_cumulative(data_->water_balance_idx, schur0->get_solution());
749  balance_->calculate_instant(data_->water_balance_idx, schur0->get_solution());
750  balance_->output();
751 }
752 
753 
755 {
756  return schur0->get_solution_precision();
757 }
758 
759 
760 // ===========================================================================================
761 //
762 // MATRIX ASSEMBLY - we use abstract assembly routine, where LS Mat/Vec SetValues
763 // are in fact pointers to allocating or filling functions - this is governed by Linsystem roitunes
764 //
765 // =======================================================================================
767 {
768  START_TIMER("DarcyFlowMH_Steady::assembly_steady_mh_matrix");
769 
770  // set auxiliary flag for switchting Dirichlet like BC
771  data_->force_no_neumann_bc = use_steady_assembly_ && (nonlinear_iteration_ == 0);
772  data_->n_schur_compls = n_schur_compls;
773 
774 
775  balance_->start_flux_assembly(data_->water_balance_idx);
776 
777  // TODO: try to move this into balance, or have it in the generic assembler class, that should perform the cell loop
778  // including various pre- and post-actions
779  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
780  unsigned int dim = dh_cell.dim();
781  assembler[dim-1]->assemble(dh_cell);
782  }
783 
784 
785  balance_->finish_flux_assembly(data_->water_balance_idx);
786 
787 }
788 
789 
791 {
792  START_TIMER("DarcyFlowMH_Steady::allocate_mh_matrix");
793 
794  // set auxiliary flag for switchting Dirichlet like BC
795  data_->n_schur_compls = n_schur_compls;
796  LinSys *ls = schur0;
797 
798  // to make space for second schur complement, max. 10 neighbour edges of one el.
799  double zeros[100000];
800  for(int i=0; i<100000; i++) zeros[i] = 0.0;
801 
802  std::vector<LongIdx> tmp_rows;
803  tmp_rows.reserve(200);
804 
805  std::vector<LongIdx> dofs, dofs_ngh;
806  dofs.reserve(data_->dh_->max_elem_dofs());
807  dofs_ngh.reserve(data_->dh_->max_elem_dofs());
808 
809  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
810  ElementAccessor<3> ele = dh_cell.elm();
811 
812  const uint ndofs = dh_cell.n_dofs();
813  dofs.resize(ndofs);
814  dh_cell.get_dof_indices(dofs);
815 
816  // whole local MH matrix
817  ls->mat_set_values(ndofs, dofs.data(), ndofs, dofs.data(), zeros);
818 
819  tmp_rows.clear();
820 
821  // compatible neighborings rows
822  unsigned int n_neighs = ele->n_neighs_vb();
823  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
824  // every compatible connection adds a 2x2 matrix involving
825  // current element pressure and a connected edge pressure
826 
827  // read neighbor dofs (dh dofhandler)
828  // neighbor cell owning neighb_side
829  DHCellAccessor dh_neighb_cell = neighb_side.cell();
830 
831  const uint ndofs_ngh = dh_neighb_cell.n_dofs();
832  dofs_ngh.resize(ndofs_ngh);
833  dh_neighb_cell.get_dof_indices(dofs_ngh);
834 
835  // local index of pedge dof on neighboring cell
836  // (dim+1) is number of edges of higher dim element
837  // TODO: replace with DHCell getter when available for FESystem component
838  const unsigned int t = dh_neighb_cell.n_dofs() - (dh_neighb_cell.dim()+1) + neighb_side.side().side_idx();
839  tmp_rows.push_back(dofs_ngh[t]);
840  }
841 
842  const uint nsides = ele->n_sides();
843  LongIdx * edge_rows = dofs.data() + nsides; // pointer to start of ele
844  // allocate always also for schur 2
845  ls->mat_set_values(nsides+1, edge_rows, n_neighs, tmp_rows.data(), zeros); // (edges, ele) x (neigh edges)
846  ls->mat_set_values(n_neighs, tmp_rows.data(), nsides+1, edge_rows, zeros); // (neigh edges) x (edges, ele)
847  ls->mat_set_values(n_neighs, tmp_rows.data(), n_neighs, tmp_rows.data(), zeros); // (neigh edges) x (neigh edges)
848 
849  tmp_rows.clear();
850 
851  if (data_->mortar_method_ != NoMortar) {
852  auto &isec_list = mesh_->mixed_intersections().element_intersections_[ele.idx()];
853  for(auto &isec : isec_list ) {
854  IntersectionLocalBase *local = isec.second;
855  DHCellAccessor dh_cell_slave = data_->dh_->cell_accessor_from_element(local->bulk_ele_idx());
856 
857  const uint ndofs_slave = dh_cell_slave.n_dofs();
858  dofs_ngh.resize(ndofs_slave);
859  dh_cell_slave.get_dof_indices(dofs_ngh);
860 
861  //DebugOut().fmt("Alloc: {} {}", ele.idx(), local->bulk_ele_idx());
862  for(unsigned int i_side=0; i_side < dh_cell_slave.elm()->n_sides(); i_side++) {
863  // TODO: replace with DHCell getter when available for FESystem component
864  tmp_rows.push_back( dofs_ngh[(ndofs_slave+1)/2+i_side] );
865  //DebugOut() << "aedge" << print_var(tmp_rows[tmp_rows.size()-1]);
866  }
867  }
868  }
869  /*
870  for(unsigned int i_side=0; i_side < ele->n_sides(); i_side++) {
871  DebugOut() << "aedge:" << print_var(edge_rows[i_side]);
872  }*/
873 
874  edge_rows = dofs.data() + nsides +1; // pointer to start of edges
875  ls->mat_set_values(nsides, edge_rows, tmp_rows.size(), tmp_rows.data(), zeros); // master edges x neigh edges
876  ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), nsides, edge_rows, zeros); // neigh edges x master edges
877  ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), tmp_rows.size(), tmp_rows.data(), zeros); // neigh edges x neigh edges
878 
879  }
880 /*
881  // alloc edge diagonal entries
882  if(rank == 0)
883  for( vector<Edge>::iterator edg = mesh_->edges.begin(); edg != mesh_->edges.end(); ++edg) {
884  int edg_idx = mh_dh.row_4_edge[edg->side(0)->edge_idx()];
885 
886 // for( vector<Edge>::iterator edg2 = mesh_->edges.begin(); edg2 != mesh_->edges.end(); ++edg2){
887 // int edg_idx2 = mh_dh.row_4_edge[edg2->side(0)->edge_idx()];
888 // if(edg_idx == edg_idx2){
889 // DBGCOUT(<< "P[ " << rank << " ] " << "edg alloc: " << edg_idx << " " << edg_idx2 << "\n");
890  ls->mat_set_value(edg_idx, edg_idx, 0.0);
891 // }
892 // }
893  }
894  */
895  /*
896  if (mortar_method_ == MortarP0) {
897  P0_CouplingAssembler(*this).assembly(*ls);
898  } else if (mortar_method_ == MortarP1) {
899  P1_CouplingAssembler(*this).assembly(*ls);
900  }*/
901 }
902 
904 {
905  START_TIMER("assembly source term");
906  balance_->start_source_assembly(data_->water_balance_idx);
907 
908  std::vector<LongIdx> global_dofs(data_->dh_->max_elem_dofs());
909 
910  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
911  ElementAccessor<3> ele = dh_cell.elm();
912 
913  const uint ndofs = dh_cell.n_dofs();
914  global_dofs.resize(ndofs);
915  dh_cell.get_dof_indices(global_dofs);
916 
917  double cs = data_->cross_section.value(ele.centre(), ele);
918 
919  // set sources
920  double source = ele.measure() * cs *
921  (data_->water_source_density.value(ele.centre(), ele)
922  +data_->extra_source.value(ele.centre(), ele));
923  // TODO: replace with DHCell getter when available for FESystem component
924  schur0->rhs_set_value(global_dofs[ndofs/2], -1.0 * source );
925 
926  balance_->add_source_values(data_->water_balance_idx, ele.region().bulk_idx(),
927  {dh_cell.get_loc_dof_indices()[ndofs/2]}, {0}, {source});
928  }
929 
930  balance_->finish_source_assembly(data_->water_balance_idx);
931 }
932 
933 
934 
935 
936 /*******************************************************************************
937  * COMPOSE WATER MH MATRIX WITHOUT SCHUR COMPLEMENT
938  ******************************************************************************/
939 
941 
942  START_TIMER("preallocation");
943 
944  if (schur0 == NULL) { // create Linear System for MH matrix
945 
946  if (in_rec.type() == LinSys_BDDC::get_input_type()) {
947 #ifdef FLOW123D_HAVE_BDDCML
948  WarningOut() << "For BDDC no Schur complements are used.";
949  n_schur_compls = 0;
950  LinSys_BDDC *ls = new LinSys_BDDC(&(*data_->dh_->distr()),
951  true); // swap signs of matrix and rhs to make the matrix SPD
952  ls->set_from_input(in_rec);
953  ls->set_solution( ele_flux_ptr->get_data_vec().petsc_vec() );
954  // possible initialization particular to BDDC
955  START_TIMER("BDDC set mesh data");
957  schur0=ls;
958  END_TIMER("BDDC set mesh data");
959 #else
960  Exception
961  xprintf(Err, "Flow123d was not build with BDDCML support.\n");
962 #endif // FLOW123D_HAVE_BDDCML
963  }
964  else if (in_rec.type() == LinSys_PETSC::get_input_type()) {
965  // use PETSC for serial case even when user wants BDDC
966  if (n_schur_compls > 2) {
967  WarningOut() << "Invalid number of Schur Complements. Using 2.";
968  n_schur_compls = 2;
969  }
970 
971  LinSys_PETSC *schur1, *schur2;
972 
973  if (n_schur_compls == 0) {
974  LinSys_PETSC *ls = new LinSys_PETSC( &(*data_->dh_->distr()) );
975 
976  // temporary solution; we have to set precision also for sequantial case of BDDC
977  // final solution should be probably call of direct solver for oneproc case
978  if (in_rec.type() != LinSys_BDDC::get_input_type()) ls->set_from_input(in_rec);
979  else {
980  ls->LinSys::set_from_input(in_rec); // get only common options
981  }
982 
983  ls->set_solution( ele_flux_ptr->get_data_vec().petsc_vec() );
984  schur0=ls;
985  } else {
986  IS is;
987  auto side_dofs_vec = get_component_indices_vec(0);
988 
989  ISCreateGeneral(PETSC_COMM_SELF, side_dofs_vec.size(), &(side_dofs_vec[0]), PETSC_COPY_VALUES, &is);
990  //ISView(is, PETSC_VIEWER_STDOUT_SELF);
991  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
992 
993  SchurComplement *ls = new SchurComplement(&(*data_->dh_->distr()), is);
994 
995  // make schur1
997  if (n_schur_compls==1) {
998  schur1 = new LinSys_PETSC(ds);
999  schur1->set_positive_definite();
1000  } else {
1001  IS is;
1002  auto elem_dofs_vec = get_component_indices_vec(1);
1003 
1004  const PetscInt *b_indices;
1005  ISGetIndices(ls->IsB, &b_indices);
1006  uint b_size = ls->loc_size_B;
1007  for(uint i_b=0, i_bb=0; i_b < b_size && i_bb < elem_dofs_vec.size(); i_b++) {
1008  if (b_indices[i_b] == elem_dofs_vec[i_bb])
1009  elem_dofs_vec[i_bb++] = i_b + ds->begin();
1010  }
1011  ISRestoreIndices(ls->IsB, &b_indices);
1012 
1013 
1014  ISCreateGeneral(PETSC_COMM_SELF, elem_dofs_vec.size(), &(elem_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1015  //ISView(is, PETSC_VIEWER_STDOUT_SELF);
1016  //OLD_ASSERT(err == 0,"Error in ISCreateStride.");
1017  SchurComplement *ls1 = new SchurComplement(ds, is); // is is deallocated by SchurComplement
1018  ls1->set_negative_definite();
1019 
1020  // make schur2
1021  schur2 = new LinSys_PETSC( ls1->make_complement_distribution() );
1022  schur2->set_positive_definite();
1023  ls1->set_complement( schur2 );
1024  schur1 = ls1;
1025  }
1026  ls->set_complement( schur1 );
1027  ls->set_from_input(in_rec);
1028  ls->set_solution( ele_flux_ptr->get_data_vec().petsc_vec() );
1029  schur0=ls;
1030  }
1031 
1032  START_TIMER("PETSC PREALLOCATION");
1033  schur0->set_symmetric();
1035 
1037 
1038  VecZeroEntries(schur0->get_solution());
1039  END_TIMER("PETSC PREALLOCATION");
1040  }
1041  else {
1042  xprintf(Err, "Unknown solver type. Internal error.\n");
1043  }
1044 
1045  END_TIMER("preallocation");
1046  }
1047 
1048 }
1049 
1050 
1051 
1052 
1054  START_TIMER("DarcyFlowMH_Steady::assembly_linear_system");
1055 
1056  data_->is_linear=true;
1057  bool is_steady = zero_time_term();
1058  //DebugOut() << "Assembly linear system\n";
1059  if (data_changed_) {
1060  data_changed_ = false;
1061  //DebugOut() << "Data changed\n";
1062  // currently we have no optimization for cases when just time term data or RHS data are changed
1063  START_TIMER("full assembly");
1064  if (typeid(*schur0) != typeid(LinSys_BDDC)) {
1065  schur0->start_add_assembly(); // finish allocation and create matrix
1066  }
1067 
1070 
1072 
1073 
1074  assembly_mh_matrix( data_->multidim_assembler ); // fill matrix
1075 
1077 // print_matlab_matrix("matrix");
1079  //MatView( *const_cast<Mat*>(schur0->get_matrix()), PETSC_VIEWER_STDOUT_WORLD );
1080  //VecView( *const_cast<Vec*>(schur0->get_rhs()), PETSC_VIEWER_STDOUT_WORLD);
1081 
1082  if (! is_steady) {
1083  START_TIMER("fix time term");
1084  //DebugOut() << "setup time term\n";
1085  // assembly time term and rhs
1086  setup_time_term();
1087  modify_system();
1088  }
1089  else
1090  {
1091  balance_->start_mass_assembly(data_->water_balance_idx);
1092  balance_->finish_mass_assembly(data_->water_balance_idx);
1093  }
1094  END_TIMER("full assembly");
1095  } else {
1096  START_TIMER("modify system");
1097  if (! is_steady) {
1098  modify_system();
1099  } else {
1100  //xprintf(PrgErr, "Planned computation time for steady solver, but data are not changed.\n");
1101  }
1102  END_TIMER("modiffy system");
1103  }
1104 
1105 }
1106 
1107 
1108 void DarcyMH::print_matlab_matrix(std::string matlab_file)
1109 {
1110  std::string output_file;
1111 
1112  if ( typeid(*schur0) == typeid(LinSys_BDDC) ){
1113 // WarningOut() << "Can output matrix only on a single processor.";
1114 // output_file = FilePath(matlab_file + "_bddc.m", FilePath::output_file);
1115 // ofstream os( output_file );
1116 // auto bddc = static_cast<LinSys_BDDC*>(schur0);
1117 // bddc->print_matrix(os);
1118  }
1119  else {//if ( typeid(*schur0) == typeid(LinSys_PETSC) ){
1120  output_file = FilePath(matlab_file + ".m", FilePath::output_file);
1121  PetscViewer viewer;
1122  PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
1123  PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1124  MatView( *const_cast<Mat*>(schur0->get_matrix()), viewer);
1125  VecView( *const_cast<Vec*>(schur0->get_rhs()), viewer);
1126  VecView( *const_cast<Vec*>(schur0->get_rhs()), viewer);
1127  VecView( *const_cast<Vec*>(&(schur0->get_solution())), viewer);
1128  }
1129 // else{
1130 // WarningOut() << "No matrix output available for the current solver.";
1131 // return;
1132 // }
1133 
1134  // compute h_min for different dimensions
1135  double d_max = std::numeric_limits<double>::max();
1136  double h1 = d_max, h2 = d_max, h3 = d_max;
1137  double he2 = d_max, he3 = d_max;
1138  for (auto ele : mesh_->elements_range()) {
1139  switch(ele->dim()){
1140  case 1: h1 = std::min(h1,ele.measure()); break;
1141  case 2: h2 = std::min(h2,ele.measure()); break;
1142  case 3: h3 = std::min(h3,ele.measure()); break;
1143  }
1144 
1145  for (unsigned int j=0; j<ele->n_sides(); j++) {
1146  switch(ele->dim()){
1147  case 2: he2 = std::min(he2, ele.side(j)->measure()); break;
1148  case 3: he3 = std::min(he3, ele.side(j)->measure()); break;
1149  }
1150  }
1151  }
1152  if(h1 == d_max) h1 = 0;
1153  if(h2 == d_max) h2 = 0;
1154  if(h3 == d_max) h3 = 0;
1155  if(he2 == d_max) he2 = 0;
1156  if(he3 == d_max) he3 = 0;
1157 
1158  FILE * file;
1159  file = fopen(output_file.c_str(),"a");
1160  fprintf(file, "nA = %d;\n", data_->dh_cr_disc_->distr()->size());
1161  fprintf(file, "nB = %d;\n", data_->dh_->mesh()->get_el_ds()->size());
1162  fprintf(file, "nBF = %d;\n", data_->dh_cr_->distr()->size());
1163  fprintf(file, "h1 = %e;\nh2 = %e;\nh3 = %e;\n", h1, h2, h3);
1164  fprintf(file, "he2 = %e;\nhe3 = %e;\n", he2, he3);
1165  fclose(file);
1166 }
1167 
1168 
1169 //template <int dim>
1170 //std::vector<arma::vec3> dof_points(DHCellAccessor cell, const Mapping<dim, 3> &mapping) {
1171 //
1172 //
1173 // vector<arma::vec::fixed<dim+1>> bary_dof_points = cell->fe()->dof_points();
1174 //
1175 // std::vector<arma::vec3> points(20);
1176 // points.resize(0);
1177 //
1178 //}
1179 //
1180 
1182  START_TIMER("DarcyFlowMH_Steady::set_mesh_data_for_bddc");
1183  // prepare mesh for BDDCML
1184  // initialize arrays
1185  // auxiliary map for creating coordinates of local dofs and global-to-local numbering
1186  std::map<int, arma::vec3> localDofMap;
1187  // connectivity for the subdomain, i.e. global dof numbers on element, stored element-by-element
1188  // Indices of Nodes on Elements
1189  std::vector<int> inet;
1190  // number of degrees of freedom on elements - determines elementwise chunks of INET array
1191  // Number of Nodes on Elements
1192  std::vector<int> nnet;
1193  // Indices of Subdomain Elements in Global Numbering - for local elements, their global indices
1194  std::vector<int> isegn;
1195 
1196  // This array is currently not used in BDDCML, it was used as an interface scaling alternative to scaling
1197  // by diagonal. It corresponds to the rho-scaling.
1198  std::vector<double> element_permeability;
1199 
1200  // maximal and minimal dimension of elements
1201  uint elDimMax = 1;
1202  uint elDimMin = 3;
1203  std::vector<LongIdx> cell_dofs_global(10);
1204 
1205 
1206 
1207  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1208  // for each element, create local numbering of dofs as fluxes (sides), pressure (element centre), Lagrange multipliers (edges), compatible connections
1209 
1210  dh_cell.get_dof_indices(cell_dofs_global);
1211 
1212  inet.insert(inet.end(), cell_dofs_global.begin(), cell_dofs_global.end());
1213  uint n_inet = cell_dofs_global.size();
1214 
1215 
1216  uint dim = dh_cell.elm().dim();
1217  elDimMax = std::max( elDimMax, dim );
1218  elDimMin = std::min( elDimMin, dim );
1219 
1220  // TODO: this is consistent with previous implementation, but may be wrong as it use global element numbering
1221  // used in sequential mesh, do global numbering of distributed elements.
1222  isegn.push_back( dh_cell.elm_idx());
1223 
1224  // TODO: use FiniteElement::dof_points
1225  for (unsigned int si=0; si<dh_cell.elm()->n_sides(); si++) {
1226  arma::vec3 coord = dh_cell.elm().side(si)->centre();
1227  // TODO: replace with DHCell getter when available for FESystem component
1228  // flux dof points
1229  localDofMap.insert( std::make_pair( cell_dofs_global[si], coord ) );
1230  // pressure trace dof points
1231  localDofMap.insert( std::make_pair( cell_dofs_global[si+dim+2], coord ) );
1232  }
1233  // pressure dof points
1234  arma::vec3 elm_centre = dh_cell.elm().centre();
1235  localDofMap.insert( std::make_pair( cell_dofs_global[dim+1], elm_centre ) );
1236 
1237  // insert dofs related to compatible connections
1238  //const Element *ele = dh_cell.elm().element();
1239  for(DHCellSide side : dh_cell.neighb_sides()) {
1240  uint neigh_dim = side.cell().elm().dim();
1241  side.cell().get_dof_indices(cell_dofs_global);
1242  int edge_row = cell_dofs_global[neigh_dim+2+side.side_idx()];
1243  localDofMap.insert( std::make_pair( edge_row, side.centre() ) );
1244  inet.push_back( edge_row );
1245  n_inet++;
1246  }
1247  nnet.push_back(n_inet);
1248 
1249 
1250  // version for rho scaling
1251  // trace computation
1252  double conduct = data_->conductivity.value( elm_centre , dh_cell.elm() );
1253  auto aniso = data_->anisotropy.value( elm_centre , dh_cell.elm() );
1254 
1255  // compute mean on the diagonal
1256  double coef = 0.;
1257  for ( int i = 0; i < 3; i++) {
1258  coef = coef + aniso.at(i,i);
1259  }
1260  // Maybe divide by cs
1261  coef = conduct*coef / 3;
1262 
1263  OLD_ASSERT( coef > 0.,
1264  "Zero coefficient of hydrodynamic resistance %f . \n ", coef );
1265  element_permeability.push_back( 1. / coef );
1266  }
1267 // uint i_inet = 0;
1268 // for(int n_dofs : nnet) {
1269 // DebugOut() << "nnet: " << n_dofs;
1270 // for(int j=0; j < n_dofs; j++, i_inet++) {
1271 // DebugOut() << "inet: " << inet[i_inet];
1272 // }
1273 // }
1274 
1275  auto distr = data_->dh_->distr();
1276 // for(auto pair : localDofMap) {
1277 // DebugOut().every_proc() << "r: " << distr->myp() << " gi: " << pair.first << "xyz: " << pair.second[0];
1278 //
1279 // }
1280 
1281 
1282  //convert set of dofs to vectors
1283  // number of nodes (= dofs) on the subdomain
1284  int numNodeSub = localDofMap.size();
1285  //ASSERT_EQ( (unsigned int)numNodeSub, data_->dh_->lsize() );
1286  // Indices of Subdomain Nodes in Global Numbering - for local nodes, their global indices
1287  std::vector<int> isngn( numNodeSub );
1288  // pseudo-coordinates of local nodes (i.e. dofs)
1289  // they need not be exact, they are used just for some geometrical considerations in BDDCML,
1290  // such as selection of corners maximizing area of a triangle, bounding boxes fro subdomains to
1291  // find candidate neighbours etc.
1292  std::vector<double> xyz( numNodeSub * 3 ) ;
1293  int ind = 0;
1294  std::map<int,arma::vec3>::iterator itB = localDofMap.begin();
1295  for ( ; itB != localDofMap.end(); ++itB ) {
1296  isngn[ind] = itB -> first;
1297 
1298  arma::vec3 coord = itB -> second;
1299  for ( int j = 0; j < 3; j++ ) {
1300  xyz[ j*numNodeSub + ind ] = coord[j];
1301  }
1302 
1303  ind++;
1304  }
1305  localDofMap.clear();
1306 
1307  // Number of Nodal Degrees of Freedom
1308  // nndf is trivially one - dofs coincide with nodes
1309  std::vector<int> nndf( numNodeSub, 1 );
1310 
1311  // prepare auxiliary map for renumbering nodes
1312  typedef std::map<int,int> Global2LocalMap_; //! type for storage of global to local map
1313  Global2LocalMap_ global2LocalNodeMap;
1314  for ( unsigned ind = 0; ind < isngn.size(); ++ind ) {
1315  global2LocalNodeMap.insert( std::make_pair( static_cast<unsigned>( isngn[ind] ), ind ) );
1316  }
1317 
1318  // renumber nodes in the inet array to locals
1319  int indInet = 0;
1320  for ( unsigned int iEle = 0; iEle < isegn.size(); iEle++ ) {
1321  int nne = nnet[ iEle ];
1322  for ( int ien = 0; ien < nne; ien++ ) {
1323 
1324  int indGlob = inet[indInet];
1325  // map it to local node
1326  Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1327  OLD_ASSERT( pos != global2LocalNodeMap.end(),
1328  "Cannot remap node index %d to local indices. \n ", indGlob );
1329  int indLoc = static_cast<int> ( pos -> second );
1330 
1331  // store the node
1332  inet[ indInet++ ] = indLoc;
1333  }
1334  }
1335 
1336  int numNodes = size;
1337  int numDofsInt = size;
1338  int spaceDim = 3; // TODO: what is the proper value here?
1339  int meshDim = elDimMax;
1340 
1341  /**
1342  * We need:
1343  * - local to global element map (possibly mesh->el_4_loc
1344  * - inet, nnet - local dof numbers per element, local numbering of only those dofs that are on owned elements
1345  * 1. collect DH local dof indices on elements, manage map from DH local indices to BDDC local dof indices
1346  * 2. map collected DH indices to BDDC indices using the map
1347  * - local BDDC dofs to global dofs, use DH to BDDC map with DH local to global map
1348  * - XYZ - permuted, collect in main loop into array of size of all DH local dofs, compress and rearrange latter
1349  * - element_permeability - in main loop
1350  */
1351  bddc_ls -> load_mesh( LinSys_BDDC::BDDCMatrixType::SPD_VIA_SYMMETRICGENERAL, spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1352 }
1353 
1354 
1355 
1356 
1357 //=============================================================================
1358 // DESTROY WATER MH SYSTEM STRUCTURE
1359 //=============================================================================
1361  if (previous_solution != nullptr) chkerr(VecDestroy(&previous_solution));
1362  if (steady_diagonal != nullptr) chkerr(VecDestroy(&steady_diagonal));
1363  if (new_diagonal != nullptr) chkerr(VecDestroy(&new_diagonal));
1364  if (steady_rhs != nullptr) chkerr(VecDestroy(&steady_rhs));
1365 
1366 
1367  if (schur0 != NULL) {
1368  delete schur0;
1369  }
1370 
1371  if (output_object) delete output_object;
1372 
1373  if(time_ != nullptr)
1374  delete time_;
1375 
1376 }
1377 
1378 
1379 /*
1380 void mat_count_off_proc_values(Mat m, Vec v) {
1381  int n, first, last;
1382  const PetscInt *cols;
1383  Distribution distr(v);
1384 
1385  int n_off = 0;
1386  int n_on = 0;
1387  int n_off_rows = 0;
1388  MatGetOwnershipRange(m, &first, &last);
1389  for (int row = first; row < last; row++) {
1390  MatGetRow(m, row, &n, &cols, PETSC_NULL);
1391  bool exists_off = false;
1392  for (int i = 0; i < n; i++)
1393  if (distr.get_proc(cols[i]) != distr.myp())
1394  n_off++, exists_off = true;
1395  else
1396  n_on++;
1397  if (exists_off)
1398  n_off_rows++;
1399  MatRestoreRow(m, row, &n, &cols, PETSC_NULL);
1400  }
1401 }
1402 */
1403 
1404 
1405 
1406 
1407 
1408 
1409 
1410 
1411 
1412 
1413 
1414 
1416 {
1417  double *local_sol = schur0->get_solution_array();
1418 
1419  // cycle over local element rows
1420 
1421  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1422  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1423  ElementAccessor<3> ele = dh_cell.elm();
1424  // set initial condition
1425  // TODO: replace with DHCell getter when available for FESystem component
1426  const IntIdx p_ele_dof = dh_cell.get_loc_dof_indices()[dh_cell.n_dofs()/2];
1427  local_sol[p_ele_dof] = data_->init_pressure.value(ele.centre(),ele);
1428  }
1429 }
1430 
1432  // save diagonal of steady matrix
1433  MatGetDiagonal(*( schur0->get_matrix() ), steady_diagonal);
1434  // save RHS
1435  VecCopy(*( schur0->get_rhs()), steady_rhs);
1436 
1437 
1438  PetscScalar *local_diagonal;
1439  VecGetArray(new_diagonal,& local_diagonal);
1440 
1441  DebugOut().fmt("Setup with dt: {}\n", time_->dt());
1442 
1443  balance_->start_mass_assembly(data_->water_balance_idx);
1444 
1445  std::vector<LongIdx> dofs;
1446  dofs.reserve(data_->dh_->max_elem_dofs());
1447 
1448  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1449  ElementAccessor<3> ele = dh_cell.elm();
1450  dofs.resize(dh_cell.n_dofs());
1451  dh_cell.get_dof_indices(dofs);
1452 
1453  // TODO: replace with DHCell getter when available for FESystem component
1454  const uint p_ele_dof = dh_cell.n_dofs() / 2;
1455  // set new diagonal
1456  double diagonal_coeff = data_->cross_section.value(ele.centre(), ele)
1457  * ( data_->storativity.value(ele.centre(), ele)
1458  +data_->extra_storativity.value(ele.centre(), ele))
1459  * ele.measure();
1460  local_diagonal[dh_cell.get_loc_dof_indices()[p_ele_dof]]= - diagonal_coeff / time_->dt();
1461 
1462  balance_->add_mass_values(data_->water_balance_idx, dh_cell,
1463  {dh_cell.get_loc_dof_indices()[p_ele_dof]},
1464  {diagonal_coeff}, 0.0);
1465  }
1466  VecRestoreArray(new_diagonal,& local_diagonal);
1467  MatDiagonalSet(*( schur0->get_matrix() ), new_diagonal, ADD_VALUES);
1468 
1470 
1471  balance_->finish_mass_assembly(data_->water_balance_idx);
1472 }
1473 
1475  START_TIMER("modify system");
1476  if (time_->is_changed_dt() && time_->step().index()>0) {
1477  double scale_factor=time_->step(-2).length()/time_->step().length();
1478  if (scale_factor != 1.0) {
1479  // if time step has changed and setup_time_term not called
1480  MatDiagonalSet(*( schur0->get_matrix() ),steady_diagonal, INSERT_VALUES);
1481 
1482  VecScale(new_diagonal, time_->last_dt()/time_->dt());
1483  MatDiagonalSet(*( schur0->get_matrix() ),new_diagonal, ADD_VALUES);
1485  }
1486  }
1487 
1488  // modify RHS - add previous solution
1489  VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, previous_solution);
1490 // VecPointwiseMult(*( schur0->get_rhs()), new_diagonal, schur0->get_solution());
1491  VecAXPY(*( schur0->get_rhs()), 1.0, steady_rhs);
1493 
1494  //VecSwap(previous_solution, schur0->get_solution());
1495 }
1496 
1497 
1498 std::shared_ptr< FieldFE<3, FieldValue<3>::VectorFixed> > DarcyMH::get_velocity_field() {
1499  return ele_flux_ptr;
1500 }
1501 
1502 
1503 /// Helper method fills range (min and max) of given component
1504 void dofs_range(unsigned int n_dofs, unsigned int &min, unsigned int &max, unsigned int component) {
1505  if (component==0) {
1506  min = 0;
1507  max = n_dofs/2;
1508  } else if (component==1) {
1509  min = n_dofs/2;
1510  max = (n_dofs+1)/2;
1511  } else {
1512  min = (n_dofs+1)/2;
1513  max = n_dofs;
1514  }
1515 }
1516 
1517 
1519  ASSERT_LT_DBG(component, 3).error("Invalid component!");
1520  unsigned int i, n_dofs, min, max;
1521  std::vector<int> dof_vec;
1522  std::vector<LongIdx> dof_indices(data_->dh_->max_elem_dofs());
1523  for ( DHCellAccessor dh_cell : data_->dh_->own_range() ) {
1524  n_dofs = dh_cell.get_dof_indices(dof_indices);
1525  dofs_range(n_dofs, min, max, component);
1526  for (i=min; i<max; ++i) dof_vec.push_back(dof_indices[i]);
1527  }
1528  return dof_vec;
1529 }
1530 
1531 
1532 //-----------------------------------------------------------------------------
1533 // vim: set cindent:
TimeGovernor & time()
Definition: equation.hh:151
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > get_velocity_field() override
Field< 3, FieldValue< 3 >::VectorFixed > field_ele_velocity
Field< 3, FieldValue< 3 >::Scalar > extra_source
Externally added storativity.
FieldSet * eq_data_
Definition: equation.hh:229
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:329
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:567
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:459
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:255
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > ele_flux_ptr
Field of flux in barycenter of every element.
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:490
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:346
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:36
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:541
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:292
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:220
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:1147
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:459
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
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:337
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:222
auto disable_where(const Field< spacedim, typename FieldValue< spacedim >::Enum > &control_field, const vector< FieldEnum > &value_list) -> Field &
Definition: field.impl.hh:180
#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:258
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:264
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:221
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.