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