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