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