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