Flow123d  3.9.1-787489b
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 =
92  it::Record("Flow_Darcy_LMH_Data",FieldCommon::field_descriptor_record_description("Flow_Darcy_LMH_Data") )
93  .copy_keys( DarcyLMH::EqFields().make_field_descriptor_type("Flow_Darcy_LMH_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 
106  it::Record ns_rec = Input::Type::Record("NonlinearSolver", "Non-linear solver settings.")
107  .declare_key("linear_solver", LinSys::get_input_type(), it::Default("{}"),
108  "Linear solver for MH problem.")
109  .declare_key("tolerance", it::Double(0.0), it::Default("1E-6"),
110  "Residual tolerance.")
111  .declare_key("min_it", it::Integer(0), it::Default("1"),
112  "Minimum number of iterations (linear solutions) to use.\nThis is usefull if the convergence criteria "
113  "does not characterize your goal well enough so it converges prematurely, possibly even without a single linear solution."
114  "If greater then 'max_it' the value is set to 'max_it'.")
115  .declare_key("max_it", it::Integer(0), it::Default("100"),
116  "Maximum number of iterations (linear solutions) of the non-linear solver.")
117  .declare_key("converge_on_stagnation", it::Bool(), it::Default("false"),
118  "If a stagnation of the nonlinear solver is detected the solver stops. "
119  "A divergence is reported by default, forcing the end of the simulation. By setting this flag to 'true', the solver "
120  "ends with convergence success on stagnation, but it reports warning about it.")
121  .close();
122 
124 
125  return it::Record("Flow_Darcy_LMH", "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 
137  .declare_key("output", DarcyFlowMHOutput::get_input_type(eq_fields, "Flow_Darcy_LMH"),
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 >("Flow_Darcy_LMH") +
154 
155 
156 
158 {
159 }
160 
161 
162 
164 : DarcyMH::EqData::EqData()
165 {
166 }
167 
168 
170 {
171  auto size = dh_p_->get_local_to_global_map().size();
172  save_local_system_.resize(size);
173  bc_fluxes_reconstruted.resize(size);
174  loc_system_.resize(size);
175  postprocess_solution_.resize(size);
176 }
177 
178 
180 {
181  std::fill(save_local_system_.begin(), save_local_system_.end(), false);
182  std::fill(bc_fluxes_reconstruted.begin(), bc_fluxes_reconstruted.end(), false);
183 }
184 
185 
186 
187 
188 
189 
190 //=============================================================================
191 // CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
192 // - do it in parallel:
193 // - initial distribution of elements, edges
194 //
195 /*! @brief CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL
196  *
197  * Parameters {Solver,NSchurs} number of performed Schur
198  * complements (0,1,2) for water flow MH-system
199  *
200  */
201 //=============================================================================
202 DarcyLMH::DarcyLMH(Mesh &mesh_in, const Input::Record in_rec, TimeGovernor *tm)
203 : DarcyFlowInterface(mesh_in, in_rec),
204  output_object(nullptr),
205  data_changed_(false),
206  read_init_cond_assembly_(nullptr),
207  mh_matrix_assembly_(nullptr),
209 {
210 
211  START_TIMER("Darcy constructor");
212  {
213  auto time_record = input_record_.val<Input::Record>("time");
214  if (tm == nullptr)
215  {
216  time_ = new TimeGovernor(time_record);
217  }
218  else
219  {
220  TimeGovernor tm_from_rec(time_record);
221  if (!tm_from_rec.is_default()) // is_default() == false when time record is present in input file
222  {
223  MessageOut() << "Duplicate key 'time', time in flow equation is already initialized from parent class!";
224  ASSERT_PERMANENT(false);
225  }
226  time_ = tm;
227  }
228  }
229 
230  eq_fields_ = make_shared<EqFields>();
231  eq_data_ = make_shared<EqData>();
232  this->eq_fieldset_ = eq_fields_.get();
233 
234  eq_fields_->set_mesh(*mesh_);
235 
236  eq_data_->is_linear=true;
237 
238  size = mesh_->n_elements() + mesh_->n_sides() + mesh_->n_edges();
239  eq_data_->mortar_method_= in_rec.val<MortarMethod>("mortar_method");
240  if (eq_data_->mortar_method_ != NoMortar) {
242  }
243 
244 
245  //side_ds->view( std::cout );
246  //el_ds->view( std::cout );
247  //edge_ds->view( std::cout );
248  //rows_ds->view( std::cout );
249 
250 }
251 
253 {
254  // DebugOut() << "t = " << time_->t() << " step_end " << time_->step().end() << "\n";
255  if(eq_data_->use_steady_assembly_)
256  {
257  // In steady case, the solution is computed with the data present at time t,
258  // and the steady state solution is valid until another change in data,
259  // which should correspond to time (t+dt).
260  // "The data change appears immediatly."
261  double next_t = time_->t() + time_->estimate_dt();
262  // DebugOut() << "STEADY next_t = " << next_t << "\n";
263  return next_t * (1 - 2*std::numeric_limits<double>::epsilon());
264  }
265  else
266  {
267  // In unsteady case, the solution is computed with the data present at time t,
268  // and the solution is valid at the time t+dt.
269  // "The data change does not appear immediatly, it is integrated over time interval dt."
270  // DebugOut() << "UNSTEADY\n";
271  return time_->t();
272  }
273 }
274 
276 //connecting data fields with mesh
277 {
278 
279  START_TIMER("data init");
280  eq_data_->mesh = mesh_;
281 
282  auto gravity_array = input_record_.val<Input::Array>("gravity");
283  std::vector<double> gvec;
284  gravity_array.copy_to(gvec);
285  gvec.push_back(0.0); // zero pressure shift
286  eq_data_->gravity_ = arma::vec(gvec);
287  eq_data_->gravity_vec_ = eq_data_->gravity_.subvec(0,2);
288 
289  FieldValue<3>::VectorFixed gvalue(eq_data_->gravity_vec_);
290  auto field_algo=std::make_shared<FieldConstant<3, FieldValue<3>::VectorFixed>>();
291  field_algo->set_value(gvalue);
292  eq_fields_->gravity_field.set(field_algo, 0.0);
293  eq_fields_->bc_gravity.set(field_algo, 0.0);
294 
295  eq_fields_->bc_pressure.add_factory(
296  std::make_shared<AddPotentialFactory<3, FieldValue<3>::Scalar> >
297  (eq_fields_->bc_gravity, eq_fields_->X(), eq_fields_->bc_piezo_head) );
298  eq_fields_->bc_switch_pressure.add_factory(
299  std::make_shared<AddPotentialFactory<3, FieldValue<3>::Scalar> >
300  (eq_fields_->bc_gravity, eq_fields_->X(), eq_fields_->bc_switch_piezo_head) );
301  eq_fields_->init_pressure.add_factory(
302  std::make_shared<AddPotentialFactory<3, FieldValue<3>::Scalar> >
303  (eq_fields_->gravity_field, eq_fields_->X(), eq_fields_->init_piezo_head) );
304 
305 
306  eq_fields_->set_input_list( this->input_record_.val<Input::Array>("input_fields"), *time_ );
307  // Check that the time step was set for the transient simulation.
308  if (! zero_time_term(true) && time_->is_default() ) {
309  //THROW(ExcAssertMsg());
310  //THROW(ExcMissingTimeGovernor() << input_record_.ei_address());
311  MessageOut() << "Missing the key 'time', obligatory for the transient problems." << endl;
312  ASSERT_PERMANENT(false);
313  }
314 
315  eq_fields_->mark_input_times(*time_);
316 }
317 
319 
320  { // init DOF handler for pressure fields
321 // std::shared_ptr< FiniteElement<0> > fe0_rt = std::make_shared<FE_RT0_disc<0>>();
322  std::shared_ptr< FiniteElement<1> > fe1_rt = std::make_shared<FE_RT0_disc<1>>();
323  std::shared_ptr< FiniteElement<2> > fe2_rt = std::make_shared<FE_RT0_disc<2>>();
324  std::shared_ptr< FiniteElement<3> > fe3_rt = std::make_shared<FE_RT0_disc<3>>();
325  std::shared_ptr< FiniteElement<0> > fe0_disc = std::make_shared<FE_P_disc<0>>(0);
326  std::shared_ptr< FiniteElement<1> > fe1_disc = std::make_shared<FE_P_disc<1>>(0);
327  std::shared_ptr< FiniteElement<2> > fe2_disc = std::make_shared<FE_P_disc<2>>(0);
328  std::shared_ptr< FiniteElement<3> > fe3_disc = std::make_shared<FE_P_disc<3>>(0);
329  std::shared_ptr< FiniteElement<0> > fe0_cr = std::make_shared<FE_CR<0>>();
330  std::shared_ptr< FiniteElement<1> > fe1_cr = std::make_shared<FE_CR<1>>();
331  std::shared_ptr< FiniteElement<2> > fe2_cr = std::make_shared<FE_CR<2>>();
332  std::shared_ptr< FiniteElement<3> > fe3_cr = std::make_shared<FE_CR<3>>();
333 // static FiniteElement<0> fe0_sys = FE_P_disc<0>(0); //TODO fix and use solution with FESystem<0>( {fe0_rt, fe0_disc, fe0_cr} )
334  FESystem<0> fe0_sys( {fe0_disc, fe0_disc, fe0_cr} );
335  FESystem<1> fe1_sys( {fe1_rt, fe1_disc, fe1_cr} );
336  FESystem<2> fe2_sys( {fe2_rt, fe2_disc, fe2_cr} );
337  FESystem<3> fe3_sys( {fe3_rt, fe3_disc, fe3_cr} );
338  MixedPtr<FESystem> fe_sys( std::make_shared<FESystem<0>>(fe0_sys), std::make_shared<FESystem<1>>(fe1_sys),
339  std::make_shared<FESystem<2>>(fe2_sys), std::make_shared<FESystem<3>>(fe3_sys) );
340  std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_sys);
341  eq_data_->dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
342  eq_data_->dh_->distribute_dofs(ds);
343  }
344 
345  init_eq_data();
347 
348  eq_fields_->add_coords_field();
349 
350  { // construct pressure, velocity and piezo head fields
351  uint rt_component = 0;
352  eq_data_->full_solution = eq_data_->dh_->create_vector();
353  auto ele_flux_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(eq_data_->dh_, &eq_data_->full_solution, rt_component);
354  eq_fields_->flux.set(ele_flux_ptr, 0.0);
355 
356  eq_fields_->field_ele_velocity.set(Model<3, FieldValue<3>::VectorFixed>::create(fn_mh_velocity(), eq_fields_->flux, eq_fields_->cross_section), 0.0);
357 
358  uint p_ele_component = 1;
359  auto ele_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(eq_data_->dh_, &eq_data_->full_solution, p_ele_component);
360  eq_fields_->field_ele_pressure.set(ele_pressure_ptr, 0.0);
361 
362  uint p_edge_component = 2;
363  auto edge_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(eq_data_->dh_, &eq_data_->full_solution, p_edge_component);
364  eq_fields_->field_edge_pressure.set(edge_pressure_ptr, 0.0);
365 
366  eq_fields_->field_ele_piezo_head.set(
367  Model<3, FieldValue<3>::Scalar>::create(fn_mh_piezohead(), eq_fields_->gravity_field, eq_fields_->X(), eq_fields_->field_ele_pressure),
368  0.0
369  );
370  }
371 
372  { // init DOF handlers represents element pressure DOFs
373  uint p_element_component = 1;
374  eq_data_->dh_p_ = std::make_shared<SubDOFHandlerMultiDim>(eq_data_->dh_,p_element_component);
375  }
376 
377  { // init DOF handlers represents edge DOFs
378  uint p_edge_component = 2;
379  eq_data_->dh_cr_ = std::make_shared<SubDOFHandlerMultiDim>(eq_data_->dh_,p_edge_component);
380  }
381 
382  { // init DOF handlers represents side DOFs
383  MixedPtr<FE_CR_disc> fe_cr_disc;
384  std::shared_ptr<DiscreteSpace> ds_cr_disc = std::make_shared<EqualOrderDiscreteSpace>( mesh_, fe_cr_disc);
385  eq_data_->dh_cr_disc_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
386  eq_data_->dh_cr_disc_->distribute_dofs(ds_cr_disc);
387  }
388 
389  eq_data_->init();
390 
391  // create solution vector for 2. Schur complement linear system
392 // p_edge_solution = new VectorMPI(eq_data_->dh_cr_->distr()->lsize());
393 // full_solution = new VectorMPI(eq_data_->dh_->distr()->lsize());
394  // this creates mpi vector from DoFHandler, including ghost values
395  eq_data_->p_edge_solution = eq_data_->dh_cr_->create_vector();
396  eq_data_->p_edge_solution_previous = eq_data_->dh_cr_->create_vector();
397  eq_data_->p_edge_solution_previous_time = eq_data_->dh_cr_->create_vector();
398 
399  // Initialize bc_switch_dirichlet to size of global boundary.
400  eq_data_->bc_switch_dirichlet.resize(mesh_->n_elements()+mesh_->bc_mesh()->n_elements(), 1);
401 
402 
403  eq_data_->nonlinear_iteration_=0;
405  .val<Input::Record>("nonlinear_solver")
406  .val<Input::AbstractRecord>("linear_solver");
407 
409 
410  // auxiliary set_time call since allocation assembly evaluates fields as well
411  //time_->step().use_fparser_ = true;
414 
415 
416  // initialization of balance object
417  balance_ = std::make_shared<Balance>("water", mesh_);
418  balance_->init_from_input(input_record_.val<Input::Record>("balance"), time());
419  eq_data_->water_balance_idx = balance_->add_quantity("water_volume");
420  balance_->allocate(eq_data_->dh_, 1);
421  balance_->units(UnitSI().m(3));
422 
423  eq_data_->balance = this->balance_; // eq_data_->balance is obsolete inherited from DarcyMH::EqData, will be remove
424  eq_data_->balance_ = this->balance_;
425 
426  this->initialize_asm();
427 }
428 
430 {
431  //eq_data_->multidim_assembler = AssemblyFlowBase::create< AssemblyLMH >(eq_fields_, eq_data_);
432 }
433 
434 //void DarcyLMH::read_initial_condition()
435 //{
436 // DebugOut().fmt("Read initial condition\n");
437 //
438 // for ( DHCellAccessor dh_cell : eq_data_->dh_->own_range() ) {
439 //
440 // LocDofVec p_indices = dh_cell.cell_with_other_dh(eq_data_->dh_p_.get()).get_loc_dof_indices();
441 // ASSERT_DBG(p_indices.n_elem == 1);
442 // LocDofVec l_indices = dh_cell.cell_with_other_dh(eq_data_->dh_cr_.get()).get_loc_dof_indices();
443 // ElementAccessor<3> ele = dh_cell.elm();
444 //
445 // // set initial condition
446 // double init_value = eq_fields_->init_pressure.value(ele.centre(),ele);
447 // unsigned int p_idx = eq_data_->dh_p_->parent_indices()[p_indices[0]];
448 // eq_data_->full_solution.set(p_idx, init_value);
449 //
450 // for (unsigned int i=0; i<ele->n_sides(); i++) {
451 // uint n_sides_of_edge = ele.side(i)->edge().n_sides();
452 // unsigned int l_idx = eq_data_->dh_cr_->parent_indices()[l_indices[i]];
453 // eq_data_->full_solution.add(l_idx, init_value/n_sides_of_edge);
454 //
455 // eq_data_->p_edge_solution.add(l_indices[i], init_value/n_sides_of_edge);
456 // }
457 // }
458 //
459 // initial_condition_postprocess();
460 //
461 // eq_data_->full_solution.ghost_to_local_begin();
462 // eq_data_->full_solution.ghost_to_local_end();
463 //
464 // eq_data_->p_edge_solution.ghost_to_local_begin();
465 // eq_data_->p_edge_solution.ghost_to_local_end();
466 // eq_data_->p_edge_solution_previous_time.copy_from(eq_data_->p_edge_solution);
467 //}
468 //
469 //void DarcyLMH::initial_condition_postprocess()
470 //{}
471 
473 {
474 
475  /* TODO:
476  * - Allow solution reconstruction (pressure and velocity) from initial condition on user request.
477  * - Steady solution as an intitial condition may be forced by setting inti_time =-1, and set data for the steady solver in that time.
478  * Solver should be able to switch from and to steady case depending on the zero time term.
479  */
480 
481  //time_->step().use_fparser_ = true;
483 
484  // zero_time_term means steady case
485  eq_data_->use_steady_assembly_ = zero_time_term();
486 
487  eq_data_->p_edge_solution.zero_entries();
488 
489  if (eq_data_->use_steady_assembly_) { // steady case
490  MessageOut() << "Flow zero time step - steady case\n";
491  //read_initial_condition(); // Possible solution guess for steady case.
492  solve_nonlinear(); // with right limit data
493  } else {
494  MessageOut() << "Flow zero time step - unsteady case\n";
495  eq_data_->time_step_ = time_->dt();
497  this->read_init_cond_asm();
498  accept_time_step(); // accept zero time step, i.e. initial condition
499 
500  // we reconstruct the initial solution here
501  // during the reconstruction assembly:
502  // - the balance objects are actually allocated
503  // - the full solution vector is computed
504  START_TIMER("DarcyFlowMH::reconstruct_solution_from_schur");
506  END_TIMER("DarcyFlowMH::reconstruct_solution_from_schur");
507  }
508  //solution_output(T,right_limit); // data for time T in any case
509  output_data();
510 }
511 
512 //=============================================================================
513 // COMPOSE and SOLVE WATER MH System possibly through Schur complements
514 //=============================================================================
516 {
517  START_TIMER("Solving MH system");
518 
519  time_->next_time();
520 
521  time_->view("DARCY"); //time governor information output
522 
523  solve_time_step();
524 
525  eq_data_->full_solution.local_to_ghost_begin();
526  eq_data_->full_solution.local_to_ghost_end();
527 }
528 
529 void DarcyLMH::solve_time_step(bool output)
530 {
531  //time_->step().use_fparser_ = true;
533  bool zero_time_term_from_left=zero_time_term();
534 
535  bool jump_time = eq_fields_->storativity.is_jump_time();
536  if (! zero_time_term_from_left) {
537  MessageOut() << "Flow time step - unsteady case\n";
538  // time term not treated as zero
539  // Unsteady solution up to the T.
540 
541  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
542  eq_data_->use_steady_assembly_ = false;
543 
544  solve_nonlinear(); // with left limit data
545  if(output)
547  if (jump_time) {
548  WarningOut() << "Output of solution discontinuous in time not supported yet.\n";
549  //solution_output(T, left_limit); // output use time T- delta*dt
550  //output_data();
551  }
552  }
553 
554  if (time_->is_end()) {
555  // output for unsteady case, end_time should not be the jump time
556  // but rether check that
557  if (! zero_time_term_from_left && ! jump_time && output)
558  output_data();
559  return;
560  }
561 
562  //time_->step().use_fparser_ = true;
564  bool zero_time_term_from_right=zero_time_term();
565  if (zero_time_term_from_right) {
566  MessageOut() << "Flow time step - steady case\n";
567  // this flag is necesssary for switching BC to avoid setting zero neumann on the whole boundary in the steady case
568  eq_data_->use_steady_assembly_ = true;
569  solve_nonlinear(); // with right limit data
570  if(output)
572 
573  } else if (! zero_time_term_from_left && jump_time) {
574  WarningOut() << "Discontinuous time term not supported yet.\n";
575  //solution_transfer(); // internally call set_time(T, left) and set_time(T,right) again
576  //solve_nonlinear(); // with right limit data
577  }
578  //solution_output(T,right_limit); // data for time T in any case
579  if (output)
580  output_data();
581 }
582 
583 bool DarcyLMH::zero_time_term(bool time_global) {
584  if (time_global) {
585  return (eq_fields_->storativity.input_list_size() == 0);
586  } else {
587  return eq_fields_->storativity.field_result(mesh_->region_db().get_region_set("BULK")) == result_zeros;
588  }
589 }
590 
591 
593 {
594 
596  double residual_norm = lin_sys_schur().compute_residual();
597  eq_data_->nonlinear_iteration_ = 0;
598  MessageOut().fmt("[nonlinear solver] norm of initial residual: {}\n", residual_norm);
599 
600  // Reduce is_linear flag.
601  int is_linear_common;
602  MPI_Allreduce(&(eq_data_->is_linear), &is_linear_common,1, MPI_INT ,MPI_MIN,PETSC_COMM_WORLD);
603 
604  Input::Record nl_solver_rec = input_record_.val<Input::Record>("nonlinear_solver");
605  this->tolerance_ = nl_solver_rec.val<double>("tolerance");
606  this->max_n_it_ = nl_solver_rec.val<unsigned int>("max_it");
607  this->min_n_it_ = nl_solver_rec.val<unsigned int>("min_it");
608  if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->max_n_it_;
609 
610  if (! is_linear_common) {
611  // set tolerances of the linear solver unless they are set by user.
612  lin_sys_schur().set_tolerances(0.1*this->tolerance_, 0.01*this->tolerance_, 100);
613  }
614  vector<double> convergence_history;
615 
616  while (eq_data_->nonlinear_iteration_ < this->min_n_it_ ||
617  (residual_norm > this->tolerance_ && eq_data_->nonlinear_iteration_ < this->max_n_it_ )) {
618  ASSERT_EQ( convergence_history.size(), eq_data_->nonlinear_iteration_ );
619  convergence_history.push_back(residual_norm);
620 
621  // print_matlab_matrix("matrix_" + std::to_string(time_->step().index()) + "_it_" + std::to_string(nonlinear_iteration_));
622  // stagnation test
623  if (convergence_history.size() >= 5 &&
624  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
625  convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
626  // stagnation
627  if (input_record_.val<Input::Record>("nonlinear_solver").val<bool>("converge_on_stagnation")) {
628  WarningOut().fmt("Accept solution on stagnation. Its: {} Residual: {}\n", eq_data_->nonlinear_iteration_, residual_norm);
629  break;
630  } else {
631  THROW(ExcSolverDiverge() << EI_Reason("Stagnation."));
632  }
633  }
634 
635  if (! is_linear_common){
636  eq_data_->p_edge_solution_previous.copy_from(eq_data_->p_edge_solution);
637  eq_data_->p_edge_solution_previous.local_to_ghost_begin();
638  eq_data_->p_edge_solution_previous.local_to_ghost_end();
639  }
640 
642  MessageOut().fmt("[schur solver] lin. it: {}, reason: {}, residual: {}\n",
643  si.n_iterations, si.converged_reason, lin_sys_schur().compute_residual());
644 
645  eq_data_->nonlinear_iteration_++;
646 
647  // hack to make BDDC work with empty compute_residual
648  if (is_linear_common){
649  // we want to print this info in linear (and steady) case
650  residual_norm = lin_sys_schur().compute_residual();
651  MessageOut().fmt("[nonlinear solver] lin. it: {}, reason: {}, residual: {}\n",
652  si.n_iterations, si.converged_reason, residual_norm);
653  break;
654  }
655  data_changed_=true; // force reassembly for non-linear case
656 
657  double alpha = 1; // how much of new solution
658  VecAXPBY(eq_data_->p_edge_solution.petsc_vec(), (1-alpha), alpha, eq_data_->p_edge_solution_previous.petsc_vec());
659 
660  //LogOut().fmt("Linear solver ended with reason: {} \n", si.converged_reason );
661  //ASSERT_PERMANENT_GE( si.converged_reason, 0).error("Linear solver failed to converge.\n");
663 
664  residual_norm = lin_sys_schur().compute_residual();
665  MessageOut().fmt("[nonlinear solver] it: {} lin. it: {}, reason: {}, residual: {}\n",
666  eq_data_->nonlinear_iteration_, si.n_iterations, si.converged_reason, residual_norm);
667  }
668 
669 // reconstruct_solution_from_schur(eq_data_->multidim_assembler);
670  START_TIMER("DarcyFlowMH::reconstruct_solution_from_schur");
672  END_TIMER("DarcyFlowMH::reconstruct_solution_from_schur");
673 
674  // adapt timestep
675  if (! this->zero_time_term()) {
676  double mult = 1.0;
677  if (eq_data_->nonlinear_iteration_ < 3) mult = 1.6;
678  if (eq_data_->nonlinear_iteration_ > 7) mult = 0.7;
679  time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
680  // int result = time_->set_upper_constraint(time_->dt() * mult, "Darcy adaptivity.");
681  //DebugOut().fmt("time adaptivity, res: {} it: {} m: {} dt: {} edt: {}\n", result, nonlinear_iteration_, mult, time_->dt(), time_->estimate_dt());
682  }
683 }
684 
685 
687 {
688  eq_data_->p_edge_solution_previous_time.copy_from(eq_data_->p_edge_solution);
689  eq_data_->p_edge_solution_previous_time.local_to_ghost_begin();
690  eq_data_->p_edge_solution_previous_time.local_to_ghost_end();
691 }
692 
693 
695  START_TIMER("Darcy output data");
696 
697  // print_matlab_matrix("matrix_" + std::to_string(time_->step().index()));
698 
699  //time_->view("DARCY"); //time governor information output
700  this->output_object->output();
701 
702 
703  START_TIMER("Darcy balance output");
704  balance_->calculate_cumulative(eq_data_->water_balance_idx, eq_data_->full_solution.petsc_vec());
705  balance_->calculate_instant(eq_data_->water_balance_idx, eq_data_->full_solution.petsc_vec());
706  balance_->output();
707 }
708 
709 
710 //double DarcyLMH::solution_precision() const
711 //{
712 // return eq_data_->lin_sys_schur->get_solution_precision();
713 //}
714 
715 
716 // ===========================================================================================
717 //
718 // MATRIX ASSEMBLY - we use abstract assembly routine, where LS Mat/Vec SetValues
719 // are in fact pointers to allocating or filling functions - this is governed by Linsystem roitunes
720 //
721 // =======================================================================================
722 //void DarcyLMH::assembly_mh_matrix(FMT_UNUSED MultidimAssembly& assembler)
723 //{
724 // START_TIMER("DarcyLMH::assembly_steady_mh_matrix");
725 //
726 // // DebugOut() << "assembly_mh_matrix \n";
727 // // set auxiliary flag for switchting Dirichlet like BC
728 // eq_data_->force_no_neumann_bc = eq_data_->use_steady_assembly_ && (eq_data_->nonlinear_iteration_ == 0);
729 //
730 // balance_->start_flux_assembly(eq_data_->water_balance_idx);
731 // balance_->start_source_assembly(eq_data_->water_balance_idx);
732 // balance_->start_mass_assembly(eq_data_->water_balance_idx);
733 //
734 // // TODO: try to move this into balance, or have it in the generic assembler class, that should perform the cell loop
735 // // including various pre- and post-actions
736 //// for ( DHCellAccessor dh_cell : eq_data_->dh_->own_range() ) {
737 //// unsigned int dim = dh_cell.dim();
738 //// assembler[dim-1]->assemble(dh_cell);
739 //// }
740 // this->mh_matrix_assembly_->assemble(eq_data_->dh_);
741 //
742 //
743 // balance_->finish_mass_assembly(eq_data_->water_balance_idx);
744 // balance_->finish_source_assembly(eq_data_->water_balance_idx);
745 // balance_->finish_flux_assembly(eq_data_->water_balance_idx);
746 //
747 //}
748 
749 
751 {
752  START_TIMER("DarcyLMH::allocate_mh_matrix");
753 
754  // to make space for second schur complement, max. 10 neighbour edges of one el.
755  double zeros[100000];
756  for(int i=0; i<100000; i++) zeros[i] = 0.0;
757 
758  std::vector<LongIdx> tmp_rows;
759  tmp_rows.reserve(200);
760 
761  std::vector<LongIdx> dofs, dofs_ngh;
762  dofs.reserve(eq_data_->dh_cr_->max_elem_dofs());
763  dofs_ngh.reserve(eq_data_->dh_cr_->max_elem_dofs());
764 
765  // DebugOut() << "Allocate new schur\n";
766  for ( DHCellAccessor dh_cell : eq_data_->dh_cr_->own_range() ) {
767  ElementAccessor<3> ele = dh_cell.elm();
768 
769  const uint ndofs = dh_cell.n_dofs();
770  dofs.resize(dh_cell.n_dofs());
771  dh_cell.get_dof_indices(dofs);
772 
773  int* dofs_ptr = dofs.data();
774  lin_sys_schur().mat_set_values(ndofs, dofs_ptr, ndofs, dofs_ptr, zeros);
775 
776  tmp_rows.clear();
777 
778  // compatible neighborings rows
779  unsigned int n_neighs = ele->n_neighs_vb();
780  for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
781  // every compatible connection adds a 2x2 matrix involving
782  // current element pressure and a connected edge pressure
783 
784  // read neighbor dofs (dh_cr dofhandler)
785  // neighbor cell owning neighb_side
786  DHCellAccessor dh_neighb_cell = neighb_side.cell();
787 
788  const uint ndofs_ngh = dh_neighb_cell.n_dofs();
789  dofs_ngh.resize(ndofs_ngh);
790  dh_neighb_cell.get_dof_indices(dofs_ngh);
791 
792  // local index of pedge dof on neighboring cell
793  tmp_rows.push_back(dofs_ngh[neighb_side.side().side_idx()]);
794  }
795 
796  lin_sys_schur().mat_set_values(ndofs, dofs_ptr, n_neighs, tmp_rows.data(), zeros); // (edges) x (neigh edges)
797  lin_sys_schur().mat_set_values(n_neighs, tmp_rows.data(), ndofs, dofs_ptr, zeros); // (neigh edges) x (edges)
798  lin_sys_schur().mat_set_values(n_neighs, tmp_rows.data(), n_neighs, tmp_rows.data(), zeros); // (neigh edges) x (neigh edges)
799 
800  tmp_rows.clear();
801  if (eq_data_->mortar_method_ != NoMortar) {
802  auto &isec_list = mesh_->mixed_intersections().element_intersections_[ele.idx()];
803  for(auto &isec : isec_list ) {
804  IntersectionLocalBase *local = isec.second;
805  DHCellAccessor dh_cell_slave = eq_data_->dh_cr_->cell_accessor_from_element(local->bulk_ele_idx());
806 
807  const uint ndofs_slave = dh_cell_slave.n_dofs();
808  dofs_ngh.resize(ndofs_slave);
809  dh_cell_slave.get_dof_indices(dofs_ngh);
810 
811  //DebugOut().fmt("Alloc: {} {}", ele.idx(), local->bulk_ele_idx());
812  for(unsigned int i_side=0; i_side < dh_cell_slave.elm()->n_sides(); i_side++) {
813  tmp_rows.push_back( dofs_ngh[i_side] );
814  //DebugOut() << "aedge" << print_var(tmp_rows[tmp_rows.size()-1]);
815  }
816  }
817  }
818 
819  lin_sys_schur().mat_set_values(ndofs, dofs_ptr, tmp_rows.size(), tmp_rows.data(), zeros); // master edges x slave edges
820  lin_sys_schur().mat_set_values(tmp_rows.size(), tmp_rows.data(), ndofs, dofs_ptr, zeros); // slave edges x master edges
821  lin_sys_schur().mat_set_values(tmp_rows.size(), tmp_rows.data(), tmp_rows.size(), tmp_rows.data(), zeros); // slave edges x slave edges
822  }
823  // DebugOut() << "end Allocate new schur\n";
824 
825  // int local_dofs[10];
826  // unsigned int nsides;
827  // for ( DHCellAccessor dh_cell : eq_data_->dh_->own_range() ) {
828  // LocalElementAccessorBase<3> ele_ac(dh_cell);
829  // nsides = ele_ac.n_sides();
830 
831  // //allocate at once matrix [sides,ele,edges]x[sides,ele,edges]
832  // loc_size = 1 + 2*nsides;
833  // unsigned int i_side = 0;
834 
835  // for (; i_side < nsides; i_side++) {
836  // local_dofs[i_side] = ele_ac.side_row(i_side);
837  // local_dofs[i_side+nsides] = ele_ac.edge_row(i_side);
838  // }
839  // local_dofs[i_side+nsides] = ele_ac.ele_row();
840  // int * edge_rows = local_dofs + nsides;
841  // //int ele_row = local_dofs[0];
842 
843  // // whole local MH matrix
844  // ls->mat_set_values(loc_size, local_dofs, loc_size, local_dofs, zeros);
845 
846 
847  // // compatible neighborings rows
848  // unsigned int n_neighs = ele_ac.element_accessor()->n_neighs_vb();
849  // unsigned int i=0;
850  // for ( DHCellSide neighb_side : dh_cell.neighb_sides() ) {
851  // //for (unsigned int i = 0; i < n_neighs; i++) {
852  // // every compatible connection adds a 2x2 matrix involving
853  // // current element pressure and a connected edge pressure
854  // Neighbour *ngh = ele_ac.element_accessor()->neigh_vb[i];
855  // DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element(neighb_side.elem_idx());
856  // LocalElementAccessorBase<3> acc_higher_dim( cell_higher_dim );
857  // for (unsigned int j = 0; j < neighb_side.element().dim()+1; j++)
858  // if (neighb_side.element()->edge_idx(j) == ngh->edge_idx()) {
859  // int neigh_edge_row = acc_higher_dim.edge_row(j);
860  // tmp_rows.push_back(neigh_edge_row);
861  // break;
862  // }
863  // //DebugOut() << "CC" << print_var(tmp_rows[i]);
864  // ++i;
865  // }
866 
867  // // allocate always also for schur 2
868  // ls->mat_set_values(nsides+1, edge_rows, n_neighs, tmp_rows.data(), zeros); // (edges, ele) x (neigh edges)
869  // ls->mat_set_values(n_neighs, tmp_rows.data(), nsides+1, edge_rows, zeros); // (neigh edges) x (edges, ele)
870  // ls->mat_set_values(n_neighs, tmp_rows.data(), n_neighs, tmp_rows.data(), zeros); // (neigh edges) x (neigh edges)
871 
872  // tmp_rows.clear();
873 
874  // if (eq_data_->mortar_method_ != NoMortar) {
875  // auto &isec_list = mesh_->mixed_intersections().element_intersections_[ele_ac.ele_global_idx()];
876  // for(auto &isec : isec_list ) {
877  // IntersectionLocalBase *local = isec.second;
878  // LocalElementAccessorBase<3> slave_acc( eq_data_->dh_->cell_accessor_from_element(local->bulk_ele_idx()) );
879  // //DebugOut().fmt("Alloc: {} {}", ele_ac.ele_global_idx(), local->bulk_ele_idx());
880  // for(unsigned int i_side=0; i_side < slave_acc.dim()+1; i_side++) {
881  // tmp_rows.push_back( slave_acc.edge_row(i_side) );
882  // //DebugOut() << "aedge" << print_var(tmp_rows[tmp_rows.size()-1]);
883  // }
884  // }
885  // }
886  // /*
887  // for(unsigned int i_side=0; i_side < ele_ac.element_accessor()->n_sides(); i_side++) {
888  // DebugOut() << "aedge:" << print_var(edge_rows[i_side]);
889  // }*/
890 
891  // ls->mat_set_values(nsides, edge_rows, tmp_rows.size(), tmp_rows.data(), zeros); // master edges x neigh edges
892  // ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), nsides, edge_rows, zeros); // neigh edges x master edges
893  // ls->mat_set_values(tmp_rows.size(), tmp_rows.data(), tmp_rows.size(), tmp_rows.data(), zeros); // neigh edges x neigh edges
894 
895  // }
896 /*
897  // alloc edge diagonal entries
898  if(rank == 0)
899  for( vector<Edge>::iterator edg = mesh_->edges.begin(); edg != mesh_->edges.end(); ++edg) {
900  int edg_idx = mh_dh.row_4_edge[edg->side(0)->edge_idx()];
901 
902 // for( vector<Edge>::iterator edg2 = mesh_->edges.begin(); edg2 != mesh_->edges.end(); ++edg2){
903 // int edg_idx2 = mh_dh.row_4_edge[edg2->side(0)->edge_idx()];
904 // if(edg_idx == edg_idx2){
905 // DBGCOUT(<< "P[ " << rank << " ] " << "edg alloc: " << edg_idx << " " << edg_idx2 << "\n");
906  ls->mat_set_value(edg_idx, edg_idx, 0.0);
907 // }
908 // }
909  }
910  */
911  /*
912  if (mortar_method_ == MortarP0) {
913  P0_CouplingAssembler(*this).assembly(*ls);
914  } else if (mortar_method_ == MortarP1) {
915  P1_CouplingAssembler(*this).assembly(*ls);
916  }*/
917 }
918 
919 
920 
921 /*******************************************************************************
922  * COMPOSE WATER MH MATRIX WITHOUT SCHUR COMPLEMENT
923  ******************************************************************************/
924 
926 
927  START_TIMER("preallocation");
928 
929  // if (schur0 == NULL) { // create Linear System for MH matrix
930 
931 // if (in_rec.type() == LinSys_BDDC::get_input_type()) {
932 // #ifdef FLOW123D_HAVE_BDDCML
933 // WarningOut() << "For BDDC no Schur complements are used.";
934 // n_schur_compls = 0;
935 // LinSys_BDDC *ls = new LinSys_BDDC(&(*eq_data_->dh_->distr()),
936 // true); // swap signs of matrix and rhs to make the matrix SPD
937 // ls->set_from_input(in_rec);
938 // ls->set_solution( eq_data_->full_solution.petsc_vec() );
939 // // possible initialization particular to BDDC
940 // START_TIMER("BDDC set mesh data");
941 // set_mesh_data_for_bddc(ls);
942 // schur0=ls;
943 // END_TIMER("BDDC set mesh data");
944 // #else
945 // Exception
946 // THROW( ExcBddcmlNotSupported() );
947 // #endif // FLOW123D_HAVE_BDDCML
948 // }
949 // else
950  if (in_rec.type() == LinSys_PETSC::get_input_type()) {
951  // use PETSC for serial case even when user wants BDDC
952 
953  eq_data_->lin_sys_schur = std::make_shared<LinSys_PETSC>( &(*eq_data_->dh_cr_->distr()) );
954  lin_sys_schur().set_from_input(in_rec);
956  lin_sys_schur().set_solution( eq_data_->p_edge_solution.petsc_vec() );
958 
959 // LinSys_PETSC *schur1, *schur2;
960 
961 // if (n_schur_compls == 0) {
962 // LinSys_PETSC *ls = new LinSys_PETSC( &(*eq_data_->dh_->distr()) );
963 
964 // // temporary solution; we have to set precision also for sequantial case of BDDC
965 // // final solution should be probably call of direct solver for oneproc case
966 // // if (in_rec.type() != LinSys_BDDC::get_input_type()) ls->set_from_input(in_rec);
967 // // else {
968 // // ls->LinSys::set_from_input(in_rec); // get only common options
969 // // }
970 // ls->set_from_input(in_rec);
971 
972 // // ls->set_solution( eq_data_->full_solution.petsc_vec() );
973 // schur0=ls;
974 // } else {
975 // IS is;
976 // auto side_dofs_vec = get_component_indices_vec(0);
977 
978 // ISCreateGeneral(PETSC_COMM_SELF, side_dofs_vec.size(), &(side_dofs_vec[0]), PETSC_COPY_VALUES, &is);
979 // //ISView(is, PETSC_VIEWER_STDOUT_SELF);
980 // //ASSERT_PERMANENT(err == 0).error("Error in ISCreateStride.");
981 
982 // SchurComplement *ls = new SchurComplement(&(*eq_data_->dh_->distr()), is);
983 
984 // // make schur1
985 // Distribution *ds = ls->make_complement_distribution();
986 // if (n_schur_compls==1) {
987 // schur1 = new LinSys_PETSC(ds);
988 // schur1->set_positive_definite();
989 // } else {
990 // IS is;
991 // auto elem_dofs_vec = get_component_indices_vec(1);
992 
993 // const PetscInt *b_indices;
994 // ISGetIndices(ls->IsB, &b_indices);
995 // uint b_size = ls->loc_size_B;
996 // for(uint i_b=0, i_bb=0; i_b < b_size && i_bb < elem_dofs_vec.size(); i_b++) {
997 // if (b_indices[i_b] == elem_dofs_vec[i_bb])
998 // elem_dofs_vec[i_bb++] = i_b + ds->begin();
999 // }
1000 // ISRestoreIndices(ls->IsB, &b_indices);
1001 
1002 
1003 // ISCreateGeneral(PETSC_COMM_SELF, elem_dofs_vec.size(), &(elem_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1004 // //ISView(is, PETSC_VIEWER_STDOUT_SELF);
1005 // //ASSERT_PERMANENT(err == 0).error("Error in ISCreateStride.");
1006 // SchurComplement *ls1 = new SchurComplement(ds, is); // is is deallocated by SchurComplement
1007 // ls1->set_negative_definite();
1008 
1009 // // make schur2
1010 // schur2 = new LinSys_PETSC( ls1->make_complement_distribution() );
1011 // schur2->set_positive_definite();
1012 // ls1->set_complement( schur2 );
1013 // schur1 = ls1;
1014 // }
1015 // ls->set_complement( schur1 );
1016 // ls->set_from_input(in_rec);
1017 // // ls->set_solution( eq_data_->full_solution.petsc_vec() );
1018 // schur0=ls;
1019  // }
1020 
1021  START_TIMER("PETSC PREALLOCATION");
1023 
1025 
1026  eq_data_->full_solution.zero_entries();
1027  eq_data_->p_edge_solution.zero_entries();
1028  END_TIMER("PETSC PREALLOCATION");
1029  }
1030  else {
1031  THROW( ExcUnknownSolver() );
1032  }
1033 
1034  END_TIMER("preallocation");
1035 }
1036 
1038 {}
1039 
1040 //void DarcyLMH::reconstruct_solution_from_schur(MultidimAssembly& assembler)
1041 //{
1042 // START_TIMER("DarcyFlowMH::reconstruct_solution_from_schur");
1043 //
1044 // eq_data_->full_solution.zero_entries();
1045 // eq_data_->p_edge_solution.local_to_ghost_begin();
1046 // eq_data_->p_edge_solution.local_to_ghost_end();
1047 //
1048 // balance_->start_flux_assembly(eq_data_->water_balance_idx);
1049 // balance_->start_source_assembly(eq_data_->water_balance_idx);
1050 // balance_->start_mass_assembly(eq_data_->water_balance_idx);
1051 //
1052 // for ( DHCellAccessor dh_cell : eq_data_->dh_->own_range() ) {
1053 // unsigned int dim = dh_cell.dim();
1054 // assembler[dim-1]->assemble_reconstruct(dh_cell);
1055 // }
1056 //
1057 // eq_data_->full_solution.local_to_ghost_begin();
1058 // eq_data_->full_solution.local_to_ghost_end();
1059 //
1060 // balance_->finish_mass_assembly(eq_data_->water_balance_idx);
1061 // balance_->finish_source_assembly(eq_data_->water_balance_idx);
1062 // balance_->finish_flux_assembly(eq_data_->water_balance_idx);
1063 //}
1064 
1066  START_TIMER("DarcyFlowMH::assembly_linear_system");
1067 // DebugOut() << "DarcyLMH::assembly_linear_system\n";
1068 
1069  eq_data_->p_edge_solution.local_to_ghost_begin();
1070  eq_data_->p_edge_solution.local_to_ghost_end();
1071 
1072  eq_data_->is_linear=true;
1073  //DebugOut() << "Assembly linear system\n";
1074 // if (data_changed_) {
1075 // data_changed_ = false;
1076  {
1077  //DebugOut() << "Data changed\n";
1078  // currently we have no optimization for cases when just time term data or RHS data are changed
1079  START_TIMER("full assembly");
1080 // if (typeid(*schur0) != typeid(LinSys_BDDC)) {
1081 // schur0->start_add_assembly(); // finish allocation and create matrix
1082 // schur_compl->start_add_assembly();
1083 // }
1084 
1086 
1089 
1090  eq_data_->time_step_ = time_->dt();
1091 
1092  START_TIMER("DarcyLMH::assembly_steady_mh_matrix");
1093  this->mh_matrix_assembly_->assemble(eq_data_->dh_);; // fill matrix
1094  END_TIMER("DarcyLMH::assembly_steady_mh_matrix");
1095 // assembly_mh_matrix( eq_data_->multidim_assembler ); // fill matrix
1096 
1099 
1100  // print_matlab_matrix("matrix");
1101  }
1102 }
1103 
1104 
1105 void DarcyLMH::print_matlab_matrix(std::string matlab_file)
1106 {
1107  std::string output_file;
1108 
1109  // compute h_min for different dimensions
1110  double d_max = std::numeric_limits<double>::max();
1111  double h1 = d_max, h2 = d_max, h3 = d_max;
1112  double he2 = d_max, he3 = d_max;
1113  for (auto ele : mesh_->elements_range()) {
1114  switch(ele->dim()){
1115  case 1: h1 = std::min(h1,ele.measure()); break;
1116  case 2: h2 = std::min(h2,ele.measure()); break;
1117  case 3: h3 = std::min(h3,ele.measure()); break;
1118  }
1119 
1120  for (unsigned int j=0; j<ele->n_sides(); j++) {
1121  switch(ele->dim()){
1122  case 2: he2 = std::min(he2, ele.side(j)->measure()); break;
1123  case 3: he3 = std::min(he3, ele.side(j)->measure()); break;
1124  }
1125  }
1126  }
1127  if(h1 == d_max) h1 = 0;
1128  if(h2 == d_max) h2 = 0;
1129  if(h3 == d_max) h3 = 0;
1130  if(he2 == d_max) he2 = 0;
1131  if(he3 == d_max) he3 = 0;
1132 
1133  FILE * file;
1134  file = fopen(output_file.c_str(),"a");
1135  fprintf(file, "nA = %d;\n", eq_data_->dh_cr_disc_->distr()->size());
1136  fprintf(file, "nB = %d;\n", eq_data_->dh_->mesh()->get_el_ds()->size());
1137  fprintf(file, "nBF = %d;\n", eq_data_->dh_cr_->distr()->size());
1138  fprintf(file, "h1 = %e;\nh2 = %e;\nh3 = %e;\n", h1, h2, h3);
1139  fprintf(file, "he2 = %e;\nhe3 = %e;\n", he2, he3);
1140  fclose(file);
1141 
1142  {
1143  output_file = FilePath(matlab_file + "_sch_new.m", FilePath::output_file);
1144  PetscViewer viewer;
1145  PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
1146  PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1147  MatView( *const_cast<Mat*>(lin_sys_schur().get_matrix()), viewer);
1148  VecView( *const_cast<Vec*>(lin_sys_schur().get_rhs()), viewer);
1149  VecView( *const_cast<Vec*>(&(lin_sys_schur().get_solution())), viewer);
1150  VecView( *const_cast<Vec*>(&(eq_data_->full_solution.petsc_vec())), viewer);
1151  }
1152 }
1153 
1154 
1155 //template <int dim>
1156 //std::vector<arma::vec3> dof_points(DHCellAccessor cell, const Mapping<dim, 3> &mapping) {
1157 //
1158 //
1159 // vector<arma::vec::fixed<dim+1>> bary_dof_points = cell->fe()->dof_points();
1160 //
1161 // std::vector<arma::vec3> points(20);
1162 // points.resize(0);
1163 //
1164 //}
1165 //
1166 
1167 // void DarcyLMH::set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls) {
1168 // START_TIMER("DarcyFlowMH_Steady::set_mesh_data_for_bddc");
1169 // // prepare mesh for BDDCML
1170 // // initialize arrays
1171 // // auxiliary map for creating coordinates of local dofs and global-to-local numbering
1172 // std::map<int, arma::vec3> localDofMap;
1173 // // connectivity for the subdomain, i.e. global dof numbers on element, stored element-by-element
1174 // // Indices of Nodes on Elements
1175 // std::vector<int> inet;
1176 // // number of degrees of freedom on elements - determines elementwise chunks of INET array
1177 // // Number of Nodes on Elements
1178 // std::vector<int> nnet;
1179 // // Indices of Subdomain Elements in Global Numbering - for local elements, their global indices
1180 // std::vector<int> isegn;
1181 //
1182 // // This array is currently not used in BDDCML, it was used as an interface scaling alternative to scaling
1183 // // by diagonal. It corresponds to the rho-scaling.
1184 // std::vector<double> element_permeability;
1185 //
1186 // // maximal and minimal dimension of elements
1187 // uint elDimMax = 1;
1188 // uint elDimMin = 3;
1189 // std::vector<LongIdx> cell_dofs_global(10);
1190 //
1191 //
1192 //
1193 // for ( DHCellAccessor dh_cell : eq_data_->dh_->own_range() ) {
1194 // // LocalElementAccessorBase<3> ele_ac(dh_cell);
1195 // // for each element, create local numbering of dofs as fluxes (sides), pressure (element centre), Lagrange multipliers (edges), compatible connections
1196 //
1197 // dh_cell.get_dof_indices(cell_dofs_global);
1198 //
1199 // inet.insert(inet.end(), cell_dofs_global.begin(), cell_dofs_global.end());
1200 // uint n_inet = cell_dofs_global.size();
1201 //
1202 //
1203 // uint dim = dh_cell.elm().dim();
1204 // elDimMax = std::max( elDimMax, dim );
1205 // elDimMin = std::min( elDimMin, dim );
1206 //
1207 // // TODO: this is consistent with previous implementation, but may be wrong as it use global element numbering
1208 // // used in sequential mesh, do global numbering of distributed elements.
1209 // isegn.push_back( dh_cell.elm_idx());
1210 //
1211 // // TODO: use FiniteElement::dof_points
1212 // for (unsigned int si=0; si<dh_cell.elm()->n_sides(); si++) {
1213 // arma::vec3 coord = dh_cell.elm().side(si)->centre();
1214 // // flux dof points
1215 // localDofMap.insert( std::make_pair( cell_dofs_global[si], coord ) );
1216 // // pressure trace dof points
1217 // localDofMap.insert( std::make_pair( cell_dofs_global[si+dim+2], coord ) );
1218 // }
1219 // // pressure dof points
1220 // arma::vec3 elm_centre = dh_cell.elm().centre();
1221 // localDofMap.insert( std::make_pair( cell_dofs_global[dim+1], elm_centre ) );
1222 //
1223 // // insert dofs related to compatible connections
1224 // //const Element *ele = dh_cell.elm().element();
1225 // for(DHCellSide side : dh_cell.neighb_sides()) {
1226 // uint neigh_dim = side.cell().elm().dim();
1227 // side.cell().get_dof_indices(cell_dofs_global);
1228 // int edge_row = cell_dofs_global[neigh_dim+2+side.side_idx()];
1229 // localDofMap.insert( std::make_pair( edge_row, side.centre() ) );
1230 // inet.push_back( edge_row );
1231 // n_inet++;
1232 // }
1233 // nnet.push_back(n_inet);
1234 //
1235 //
1236 // // version for rho scaling
1237 // // trace computation
1238 // double conduct = eq_fields_->conductivity.value( elm_centre , dh_cell.elm() );
1239 // auto aniso = eq_fields_->anisotropy.value( elm_centre , dh_cell.elm() );
1240 //
1241 // // compute mean on the diagonal
1242 // double coef = 0.;
1243 // for ( int i = 0; i < 3; i++) {
1244 // coef = coef + aniso.at(i,i);
1245 // }
1246 // // Maybe divide by cs
1247 // coef = conduct*coef / 3;
1248 //
1249 // ASSERT_PERMANENT_GT(coef, 0).error("Zero coefficient of hydrodynamic resistance.\n");
1250 // element_permeability.push_back( 1. / coef );
1251 // }
1252 // // uint i_inet = 0;
1253 // // for(int n_dofs : nnet) {
1254 // // DebugOut() << "nnet: " << n_dofs;
1255 // // for(int j=0; j < n_dofs; j++, i_inet++) {
1256 // // DebugOut() << "inet: " << inet[i_inet];
1257 // // }
1258 // // }
1259 //
1260 // auto distr = eq_data_->dh_->distr();
1261 // // for(auto pair : localDofMap) {
1262 // // DebugOut().every_proc() << "r: " << distr->myp() << " gi: " << pair.first << "xyz: " << pair.second[0];
1263 // //
1264 // // }
1265 //
1266 //
1267 // //convert set of dofs to vectors
1268 // // number of nodes (= dofs) on the subdomain
1269 // int numNodeSub = localDofMap.size();
1270 // //ASSERT_PERMANENT_EQ( (unsigned int)numNodeSub, eq_data_->dh_->lsize() );
1271 // // Indices of Subdomain Nodes in Global Numbering - for local nodes, their global indices
1272 // std::vector<int> isngn( numNodeSub );
1273 // // pseudo-coordinates of local nodes (i.e. dofs)
1274 // // they need not be exact, they are used just for some geometrical considerations in BDDCML,
1275 // // such as selection of corners maximizing area of a triangle, bounding boxes fro subdomains to
1276 // // find candidate neighbours etc.
1277 // std::vector<double> xyz( numNodeSub * 3 ) ;
1278 // int ind = 0;
1279 // std::map<int,arma::vec3>::iterator itB = localDofMap.begin();
1280 // for ( ; itB != localDofMap.end(); ++itB ) {
1281 // isngn[ind] = itB -> first;
1282 //
1283 // arma::vec3 coord = itB -> second;
1284 // for ( int j = 0; j < 3; j++ ) {
1285 // xyz[ j*numNodeSub + ind ] = coord[j];
1286 // }
1287 //
1288 // ind++;
1289 // }
1290 // localDofMap.clear();
1291 //
1292 // // Number of Nodal Degrees of Freedom
1293 // // nndf is trivially one - dofs coincide with nodes
1294 // std::vector<int> nndf( numNodeSub, 1 );
1295 //
1296 // // prepare auxiliary map for renumbering nodes
1297 // typedef std::map<int,int> Global2LocalMap_; //! type for storage of global to local map
1298 // Global2LocalMap_ global2LocalNodeMap;
1299 // for ( unsigned ind = 0; ind < isngn.size(); ++ind ) {
1300 // global2LocalNodeMap.insert( std::make_pair( static_cast<unsigned>( isngn[ind] ), ind ) );
1301 // }
1302 //
1303 // // renumber nodes in the inet array to locals
1304 // int indInet = 0;
1305 // for ( unsigned int iEle = 0; iEle < isegn.size(); iEle++ ) {
1306 // int nne = nnet[ iEle ];
1307 // for ( int ien = 0; ien < nne; ien++ ) {
1308 //
1309 // int indGlob = inet[indInet];
1310 // // map it to local node
1311 // Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1312 // ASSERT_PERMANENT( pos != global2LocalNodeMap.end())(indGlob).error("Cannot remap node index to local indices. \n " );
1313 // int indLoc = static_cast<int> ( pos -> second );
1314 //
1315 // // store the node
1316 // inet[ indInet++ ] = indLoc;
1317 // }
1318 // }
1319 //
1320 // int numNodes = size;
1321 // int numDofsInt = size;
1322 // int spaceDim = 3; // TODO: what is the proper value here?
1323 // int meshDim = elDimMax;
1324 //
1325 // /**
1326 // * We need:
1327 // * - local to global element map (possibly mesh->el_4_loc
1328 // * - inet, nnet - local dof numbers per element, local numbering of only those dofs that are on owned elements
1329 // * 1. collect DH local dof indices on elements, manage map from DH local indices to BDDC local dof indices
1330 // * 2. map collected DH indices to BDDC indices using the map
1331 // * - local BDDC dofs to global dofs, use DH to BDDC map with DH local to global map
1332 // * - XYZ - permuted, collect in main loop into array of size of all DH local dofs, compress and rearrange latter
1333 // * - element_permeability - in main loop
1334 // */
1335 // bddc_ls -> load_mesh( LinSys_BDDC::BDDCMatrixType::SPD_VIA_SYMMETRICGENERAL, spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1336 // }
1337 
1338 
1339 
1340 
1341 //=============================================================================
1342 // DESTROY WATER MH SYSTEM STRUCTURE
1343 //=============================================================================
1345  if (output_object) delete output_object;
1346 
1347  if(time_ != nullptr)
1348  delete time_;
1349 
1350  if (read_init_cond_assembly_!=nullptr) {
1351  delete read_init_cond_assembly_;
1352  read_init_cond_assembly_ = nullptr;
1353  }
1354  if (mh_matrix_assembly_!=nullptr) {
1355  delete mh_matrix_assembly_;
1356  mh_matrix_assembly_ = nullptr;
1357  }
1358  if (reconstruct_schur_assembly_!=nullptr) {
1360  reconstruct_schur_assembly_ = nullptr;
1361  }
1362 }
1363 
1364 
1366  ASSERT_LT(component, 3).error("Invalid component!");
1367  unsigned int i, n_dofs, min, max;
1368  std::vector<int> dof_vec;
1369  std::vector<LongIdx> dof_indices(eq_data_->dh_->max_elem_dofs());
1370  for ( DHCellAccessor dh_cell : eq_data_->dh_->own_range() ) {
1371  n_dofs = dh_cell.get_dof_indices(dof_indices);
1372  dofs_range(n_dofs, min, max, component);
1373  for (i=min; i<max; ++i) dof_vec.push_back(dof_indices[i]);
1374  }
1375  return dof_vec;
1376 }
1377 
1378 
1383 }
1384 
1385 
1387  this->read_init_cond_assembly_->assemble(eq_data_->dh_cr_);
1388 }
1389 
1390 
1391 
1392 //-----------------------------------------------------------------------------
1393 // vim: set cindent:
LimitSide::right
@ right
Element::n_neighs_vb
unsigned int n_neighs_vb() const
Return number of neighbours.
Definition: elements.h:65
EquationBase::mesh_
Mesh * mesh_
Definition: equation.hh:219
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::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:1105
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:583
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
DarcyFlowMHOutput::output
void output()
Calculate values for output.
Definition: darcy_flow_mh_output.cc:246
TimeGovernor::dt
double dt() const
Definition: time_governor.hh:565
DarcyLMH::EqData::reset
void reset()
Reset data members.
Definition: darcy_flow_lmh.cc:179
factory.hh
field_constant.hh
vector_mpi.hh
DarcyLMH::output_object
DarcyFlowMHOutput * output_object
Definition: darcy_flow_lmh.hh:334
time_governor.hh
Basic time management class.
field_add_potential.hh
DarcyLMH::EqData
Definition: darcy_flow_lmh.hh:158
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:686
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:925
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:220
assembly_lmh.hh
distribution.hh
Support classes for parallel programing.
bc_mesh.hh
DarcyLMH::DarcyFlowMHOutput
friend class DarcyFlowMHOutput
Definition: darcy_flow_lmh.hh:349
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:275
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:148
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()
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:343
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:338
system.hh
fn_mh_velocity
Definition: assembly_models.hh:29
assembly_models.hh
Functors of FieldModels used in Darcy flow module.
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:529
DarcyLMH::update_solution
void update_solution() override
Definition: darcy_flow_lmh.cc:515
DarcyLMH::EqFields
Definition: darcy_flow_lmh.hh:151
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:1379
darcy_flow_lmh.hh
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
DarcyLMH::postprocess
virtual void postprocess()
Definition: darcy_flow_lmh.cc:1037
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:1608
DarcyLMH::assembly_linear_system
virtual void assembly_linear_system()
Definition: darcy_flow_lmh.cc:1065
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:355
DarcyLMH::min_n_it_
unsigned int min_n_it_
Definition: darcy_flow_lmh.hh:342
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:204
ASSERT_LT
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
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:472
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:62
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
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
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
DarcyLMH::eq_fields
EqFields & eq_fields()
Definition: darcy_flow_lmh.hh:236
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
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:592
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:332
DarcyLMH::registrar
static const int registrar
Registrar of class to factory.
Definition: darcy_flow_lmh.hh:359
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:163
TimeGovernor::step
const TimeStep & step(int index=-1) const
Definition: time_governor.cc:756
DarcyLMH::solved_time
virtual double solved_time() override
Definition: darcy_flow_lmh.cc:252
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
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:221
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
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
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:429
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:346
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
DarcyLMH::~DarcyLMH
virtual ~DarcyLMH() override
Definition: darcy_flow_lmh.cc:1344
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
MixedMeshIntersections::element_intersections_
std::vector< std::vector< ILpair > > element_intersections_
Definition: mixed_mesh_intersections.hh:83
partitioning.hh
DarcyLMH::initialize
void initialize() override
Definition: darcy_flow_lmh.cc:318
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:1365
DarcyLMH::EqData::init
void init()
Initialize vectors, ...
Definition: darcy_flow_lmh.cc:169
DarcyLMH::lin_sys_schur
LinSys & lin_sys_schur()
Getter for the linear system of the 2. Schur complement.
Definition: darcy_flow_lmh.hh:323
LinSys::set_from_input
virtual void set_from_input(const Input::Record in_rec)
Definition: linsys.hh:640
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:341
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::allocate_mh_matrix
void allocate_mh_matrix()
Definition: darcy_flow_lmh.cc:750
EquationBase::record_template
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Definition: equation.cc:35
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::reconstruct_schur_assembly_
GenericAssemblyBase * reconstruct_schur_assembly_
Definition: darcy_flow_lmh.hh:356
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
DarcyLMH::size
int size
Definition: darcy_flow_lmh.hh:336
ElementAccessor::idx
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Definition: accessors.hh:215
FilePath::output_file
@ output_file
Definition: file_path.hh:69
EquationBase::eq_fieldset_
FieldSet * eq_fieldset_
Definition: equation.hh:228
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:354
DarcyFlowInterface::MortarP0
@ MortarP0
Definition: darcy_flow_interface.hh:32
balance.hh
local_to_global_map.hh
MeshBase::n_elements
unsigned int n_elements() const
Definition: mesh.h:111
DarcyLMH::read_init_cond_asm
virtual void read_init_cond_asm()
Call assemble of read_init_cond_assembly_.
Definition: darcy_flow_lmh.cc:1386
IntersectionLocalBase::bulk_ele_idx
unsigned int bulk_ele_idx() const
Returns index of bulk element.
Definition: intersection_local.hh:77
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
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:345
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
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
Definition: darcy_flow_mh.hh:125
DarcyLMH::output_data
virtual void output_data() override
Write computed fields.
Definition: darcy_flow_lmh.cc:694
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:202
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
END_TIMER
#define END_TIMER(tag)
Ends a timer with specified tag.
Definition: sys_profiler.hh:149
FieldValue
Definition: field_values.hh:645
range_wrapper.hh
Implementation of range helper class.
IntersectionLocalBase
Common base for intersection object.
Definition: intersection_local.hh:48
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