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