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