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