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