Flow123d  JB_transport-d4c8564
transport_dg.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 transport_dg.cc
15 * @brief Discontinuous Galerkin method for equation of transport with dispersion.
16 * @author Jan Stebel
17 */
18 
19 #include "system/index_types.hh"
20 #include "system/sys_profiler.hh"
22 
23 #include "io/output_time.hh"
25 #include "fem/fe_values.hh"
26 #include "fem/fe_p.hh"
27 #include "fem/fe_rt.hh"
28 #include "fem/dh_cell_accessor.hh"
29 #include "fields/field_fe.hh"
31 #include "la/linsys_PETSC.hh"
32 #include "coupling/balance.hh"
37 #include "transport/heat_model.hh"
38 #include "transport/assembly_dg.hh"
39 
40 #include "fields/multi_field.hh"
41 #include "fields/generic_field.hh"
42 #include "input/factory.hh"
44 #include "mesh/accessors.hh"
45 
46 FLOW123D_FORCE_LINK_IN_CHILD(concentrationTransportModel)
48 
49 
50 
51 using namespace Input::Type;
52 
53 template<class Model>
55  return Selection("DG_variant", "Type of penalty term.")
56  .add_value(non_symmetric, "non-symmetric", "non-symmetric weighted interior penalty DG method")
57  .add_value(incomplete, "incomplete", "incomplete weighted interior penalty DG method")
58  .add_value(symmetric, "symmetric", "symmetric weighted interior penalty DG method")
59  .close();
60 }
61 
62 /*
63 * Should be removed
64 template<class Model>
65 const Selection & TransportDG<Model>::EqData::get_output_selection() {
66  return Model::ModelEqData::get_output_selection_input_type(
67  "DG",
68  "Implicit in time Discontinuous Galerkin solver")
69  .copy_values(EqData().make_output_field_selection("").close())
70  ConvectionTransport::EqData().output_fields
71  .make_output_field_selection(
72  "ConvectionTransport_output_fields",
73  "Selection of output fields for Convection Solute Transport model.")
74  .close()),
75  .close();
76 }
77 */
78 
79 template<class Model>
81  std::string equation_name = std::string(Model::ModelEqData::name()) + "_DG";
82  return Model::get_input_type("DG", "Discontinuous Galerkin (DG) solver")
85  "Solver for the linear system.")
86  .declare_key("input_fields", Array(
88  .make_field_descriptor_type(equation_name)),
90  "Input fields of the equation.")
92  "Variant of the interior penalty discontinuous Galerkin method.")
93  .declare_key("dg_order", Integer(0,3), Default("1"),
94  "Polynomial order for the finite element in DG method (order 0 is suitable if there is no diffusion/dispersion).")
95  .declare_key("init_projection", Bool(), Default("true"),
96  "If true, use DG projection of the initial condition field."
97  "Otherwise, evaluate initial condition field directly (well suited for reading native data).")
98  .declare_key("output",
99  EqFields().output_fields.make_output_type(equation_name, ""),
100  IT::Default("{ \"fields\": [ " + Model::ModelEqData::default_output_field() + "] }"),
101  "Specification of output fields and output times.")
102  .close();
103 }
104 
105 template<class Model>
107  Input::register_class< TransportDG<Model>, Mesh &, const Input::Record>(std::string(Model::ModelEqData::name()) + "_DG") +
109 
110 
111 
112 template<class Model>
114 {
115  *this+=fracture_sigma
116  .name("fracture_sigma")
117  .description(
118  "Coefficient of diffusive transfer through fractures (for each substance).")
120  .input_default("1.0")
122 
123  *this+=dg_penalty
124  .name("dg_penalty")
125  .description(
126  "Penalty parameter influencing the discontinuity of the solution (for each substance). "
127  "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
129  .input_default("1.0")
131 
132  *this += region_id.name("region_id")
135  .description("Region ids.");
136 
137  *this += subdomain.name("subdomain")
140  .description("Subdomain ids of the domain decomposition.");
141 
142 
143  // add all input fields to the output list
144  output_fields += *this;
145 
146  this->add_coords_field();
147  this->set_default_fieldset();
148 }
149 
150 
151 
152 // return the ratio of longest and shortest edge
153 template<class Model>
155 {
156  double h_max = 0, h_min = numeric_limits<double>::infinity();
157  for (unsigned int i=0; i<e->n_nodes(); i++)
158  for (unsigned int j=i+1; j<e->n_nodes(); j++)
159  {
160  double dist = arma::norm(*e.node(i) - *e.node(j));
161  h_max = max(h_max, dist);
162  h_min = min(h_min, dist);
163  }
164  return h_max/h_min;
165 }
166 
167 
168 
169 template<class Model>
171  const int K_size,
172  const vector<arma::mat33> &K,
173  const double flux,
174  const arma::vec3 &normal_vector,
175  const double alpha,
176  double &gamma)
177 {
178  double delta = 0, h = 0;
179 
180  // calculate the side diameter
181  if (side.dim() == 0)
182  {
183  h = 1;
184  }
185  else
186  {
187  for (unsigned int i=0; i<side.n_nodes(); i++)
188  for (unsigned int j=i+1; j<side.n_nodes(); j++) {
189  double dist = arma::norm(*side.node(i) - *side.node(j));
190  h = max(h, dist);
191  }
192 
193  }
194 
195  // delta is set to the average value of Kn.n on the side
196  for (int k=0; k<K_size; k++)
197  delta += dot(K[k]*normal_vector,normal_vector);
198  delta /= K_size;
199 
200  gamma = 0.5*fabs(flux) + alpha/h*delta*elem_anisotropy(side.element());
201 }
202 
203 
204 
205 template<typename Model>
207  : Model(init_mesh, in_rec),
208  input_rec(in_rec),
209  allocation_done(false),
210  mass_assembly_(nullptr)
211 {
212  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
213  // due to constexpr optimization.
214  START_TIMER(Model::ModelEqData::name());
215  // Check that Model is derived from AdvectionDiffusionModel.
217 
218  eq_data_ = make_shared<EqData>();
219  eq_fields_ = make_shared<EqFields>();
220  this->eq_fieldset_ = eq_fields_;
221  Model::init_balance(in_rec);
222 
223 
224  // Set up physical parameters.
225  eq_fields_->set_mesh(init_mesh);
226  eq_fields_->region_id = GenericField<3>::region_id(*Model::mesh_);
227  eq_fields_->subdomain = GenericField<3>::subdomain(*Model::mesh_);
228 
229 
230  // DG data parameters
231  eq_data_->dg_variant = in_rec.val<DGVariant>("dg_variant");
232  eq_data_->dg_order = in_rec.val<unsigned int>("dg_order");
233 
234  Model::init_from_input(in_rec);
235 
236  MixedPtr<FE_P_disc> fe(eq_data_->dg_order);
237  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>(Model::mesh_, fe);
238  eq_data_->dh_ = make_shared<DOFHandlerMultiDim>(*Model::mesh_);
239  eq_data_->dh_->distribute_dofs(ds);
240  //DebugOut().fmt("TDG: solution size {}\n", eq_data_->dh_->n_global_dofs());
241 
242 }
243 
244 
245 template<class Model>
247 {
248  eq_fields_->set_components(eq_data_->substances_.names());
249  eq_fields_->set_input_list( input_rec.val<Input::Array>("input_fields"), *(Model::time_) );
250  eq_data_->set_time_governor(Model::time_);
251  eq_data_->balance_ = this->balance();
252  eq_fields_->initialize(input_rec);
253 
254  // DG stabilization parameters on boundary edges
255  eq_data_->gamma.resize(eq_data_->n_substances());
256  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
257  eq_data_->gamma[sbi].resize(Model::mesh_->boundary_.size());
258 
259  // Resize coefficient arrays
260  eq_data_->max_edg_sides = max(Model::mesh_->max_edge_sides(1), max(Model::mesh_->max_edge_sides(2), Model::mesh_->max_edge_sides(3)));
261  ret_sources.resize(eq_data_->n_substances());
262  ret_sources_prev.resize(eq_data_->n_substances());
263 
264  Input::Array user_fields_arr;
265  if (input_rec.opt_val("user_fields", user_fields_arr)) {
266  this->init_user_fields(user_fields_arr, eq_fields_->output_fields);
267  }
268 
269  eq_data_->output_vec.resize(eq_data_->n_substances());
270  eq_fields_->output_field.set_components(eq_data_->substances_.names());
271  eq_fields_->output_field.set_mesh(*Model::mesh_);
272  eq_fields_->output_fields.set_mesh(*Model::mesh_);
273  eq_fields_->output_type(OutputTime::CORNER_DATA);
274 
275  eq_fields_->output_field.setup_components();
276  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
277  {
278  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
279  auto output_field_ptr = create_field_fe< 3, FieldValue<3>::Scalar >(eq_data_->dh_);
280  eq_fields_->output_field[sbi].set(output_field_ptr, 0);
281  eq_data_->output_vec[sbi] = output_field_ptr->vec();
282  }
283 
284  // set time marks for writing the output
285  eq_fields_->output_fields.initialize(Model::output_stream_, Model::mesh_, input_rec.val<Input::Record>("output"), this->time());
286 
287  // equation default PETSc solver options
288  std::string petsc_default_opts;
289  if (eq_data_->dh_->distr()->np() == 1)
290  petsc_default_opts = "-ksp_type bcgs -pc_type ilu -pc_factor_levels 2 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
291  else
292  petsc_default_opts = "-ksp_type bcgs -ksp_diagonal_scale_fix -pc_type asm -pc_asm_overlap 4 -sub_pc_type ilu -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
293 
294  // allocate matrix and vector structures
295  eq_data_->ls = new LinSys*[eq_data_->n_substances()];
296  eq_data_->ls_dt = new LinSys*[eq_data_->n_substances()];
297  eq_data_->conc_fe.resize(eq_data_->n_substances());
298 
299  MixedPtr<FE_P_disc> fe(0);
300  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>(Model::mesh_, fe);
301  eq_data_->dh_p0 = make_shared<DOFHandlerMultiDim>(*Model::mesh_);
302  eq_data_->dh_p0->distribute_dofs(ds);
303 
304  stiffness_matrix.resize(eq_data_->n_substances(), nullptr);
305  mass_matrix.resize(eq_data_->n_substances(), nullptr);
306  rhs.resize(eq_data_->n_substances(), nullptr);
307  mass_vec.resize(eq_data_->n_substances(), nullptr);
308  eq_data_->ret_vec.resize(eq_data_->n_substances(), nullptr);
309 
310  for (unsigned int sbi = 0; sbi < eq_data_->n_substances(); sbi++) {
311  eq_data_->ls[sbi] = new LinSys_PETSC(eq_data_->dh_->distr().get(), petsc_default_opts);
312  ( (LinSys_PETSC *)eq_data_->ls[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
313  eq_data_->ls[sbi]->set_solution(eq_data_->output_vec[sbi].petsc_vec());
314 
315  eq_data_->ls_dt[sbi] = new LinSys_PETSC(eq_data_->dh_->distr().get(), petsc_default_opts);
316  ( (LinSys_PETSC *)eq_data_->ls_dt[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
317 
318  eq_data_->conc_fe[sbi] = create_field_fe< 3, FieldValue<3>::Scalar >(eq_data_->dh_p0);
319 
320  VecDuplicate(eq_data_->ls[sbi]->get_solution(), &eq_data_->ret_vec[sbi]);
321  }
322 
323 
324  init_projection = input_rec.val<bool>("init_projection");
325 
326  // create assemblation object, finite element structures and distribute DOFs
327  mass_assembly_ = new GenericAssembly< MassAssemblyDim >(eq_fields_.get(), eq_data_.get());
328  stiffness_assembly_ = new GenericAssembly< StiffnessAssemblyDim >(eq_fields_.get(), eq_data_.get());
329  sources_assembly_ = new GenericAssembly< SourcesAssemblyDim >(eq_fields_.get(), eq_data_.get());
330  bdr_cond_assembly_ = new GenericAssembly< BdrConditionAssemblyDim >(eq_fields_.get(), eq_data_.get());
331 
332  if(init_projection)
333  init_assembly_ = new GenericAssembly< InitProjectionAssemblyDim >(eq_fields_.get(), eq_data_.get());
334  else
335  init_assembly_ = new GenericAssembly< InitConditionAssemblyDim >(eq_fields_.get(), eq_data_.get());
336 
337  // initialization of balance object
338  Model::balance_->allocate(eq_data_->dh_->distr()->lsize(), mass_assembly_->eval_points()->max_size());
339 
340  int qsize = mass_assembly_->eval_points()->max_size();
341  eq_data_->dif_coef.resize(eq_data_->n_substances());
342  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
343  {
344  eq_data_->dif_coef[sbi].resize(qsize);
345  }
346 
347  eq_fields_->init_condition.setup_components();
348  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
349  {
350  eq_fields_->init_condition[sbi].add_factory( std::make_shared<FieldFE<3, FieldValue<3>::Scalar>::NativeFactory>(sbi, eq_data_->dh_));
351  }
352 }
353 
354 
355 template<class Model>
357 {
358  delete Model::time_;
359 
360  if (eq_data_->gamma.size() > 0) {
361  // initialize called
362 
363  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
364  {
365  if (eq_data_->ls != nullptr) {
366  delete eq_data_->ls[i];
367  delete eq_data_->ls_dt[i];
368  }
369 
370  if (stiffness_matrix.size() > 0) {
371  if (stiffness_matrix[i])
372  chkerr(MatDestroy(&stiffness_matrix[i]));
373  if (mass_matrix[i])
374  chkerr(MatDestroy(&mass_matrix[i]));
375  if (rhs[i])
376  chkerr(VecDestroy(&rhs[i]));
377  if (mass_vec[i])
378  chkerr(VecDestroy(&mass_vec[i]));
379  if (eq_data_->ret_vec[i])
380  chkerr(VecDestroy(&eq_data_->ret_vec[i]));
381  }
382  }
383  if (eq_data_->ls != nullptr) {
384  delete[] eq_data_->ls;
385  delete[] eq_data_->ls_dt;
386  eq_data_->ls = nullptr;
387  }
388  //delete[] stiffness_matrix;
389  //delete[] mass_matrix;
390  //delete[] rhs;
391  //delete[] mass_vec;
392  //delete[] ret_vec;
393 
394  if (mass_assembly_ != nullptr) {
395  delete mass_assembly_;
396  delete stiffness_assembly_;
397  delete sources_assembly_;
398  delete bdr_cond_assembly_;
399  delete init_assembly_;
400  }
401  }
402 
403 
404 }
405 
406 
407 template<class Model>
409 {
410  START_TIMER(Model::ModelEqData::name());
411  eq_fields_->mark_input_times( *(Model::time_) );
412  eq_fields_->set_time(Model::time_->step(), LimitSide::left);
413  std::stringstream ss; // print warning message with table of uninitialized fields
414  if ( FieldCommon::print_message_table(ss, "transport DG") ) {
415  WarningOut() << ss.str();
416  }
417 
418 
419  // set initial conditions
420  set_initial_condition();
421  for (unsigned int sbi = 0; sbi < eq_data_->n_substances(); sbi++)
422  ( (LinSys_PETSC *)eq_data_->ls[sbi] )->set_initial_guess_nonzero();
423 
424  // check first time assembly - needs preallocation
425  if (!allocation_done) preallocate();
426 
427  // after preallocation we assemble the matrices and vectors required for mass balance
428  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
429  {
430  Model::balance_->calculate_instant(eq_data_->subst_idx_[sbi], eq_data_->ls[sbi]->get_solution());
431 
432  // add sources due to sorption
433  ret_sources_prev[sbi] = 0;
434  }
435 
436  output_data();
437 }
438 
439 
440 template<class Model>
442 {
443  // preallocate system matrix
444  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
445  {
446  // preallocate system matrix
447  eq_data_->ls[i]->start_allocation();
448  stiffness_matrix[i] = NULL;
449  rhs[i] = NULL;
450 
451  // preallocate mass matrix
452  eq_data_->ls_dt[i]->start_allocation();
453  mass_matrix[i] = NULL;
454  VecZeroEntries(eq_data_->ret_vec[i]);
455  }
456  stiffness_assembly_->assemble(eq_data_->dh_);
457  mass_assembly_->assemble(eq_data_->dh_);
458  sources_assembly_->assemble(eq_data_->dh_);
459  bdr_cond_assembly_->assemble(eq_data_->dh_);
460  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
461  {
462  VecAssemblyBegin(eq_data_->ret_vec[i]);
463  VecAssemblyEnd(eq_data_->ret_vec[i]);
464  }
465 
466  allocation_done = true;
467 }
468 
469 
470 
471 template<class Model>
473 {
474  START_TIMER("DG-ONE STEP");
475 
476  Model::time_->next_time();
477  Model::time_->view("TDG");
478 
479  START_TIMER("data reinit");
480  eq_fields_->set_time(Model::time_->step(), LimitSide::left);
481  END_TIMER("data reinit");
482 
483  // assemble mass matrix
484  if (mass_matrix[0] == NULL || eq_fields_->subset(FieldFlag::in_time_term).changed() )
485  {
486  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
487  {
488  eq_data_->ls_dt[i]->start_add_assembly();
489  eq_data_->ls_dt[i]->mat_zero_entries();
490  VecZeroEntries(eq_data_->ret_vec[i]);
491  }
492  mass_assembly_->assemble(eq_data_->dh_);
493  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
494  {
495  eq_data_->ls_dt[i]->finish_assembly();
496  VecAssemblyBegin(eq_data_->ret_vec[i]);
497  VecAssemblyEnd(eq_data_->ret_vec[i]);
498  // construct mass_vec for initial time
499  if (mass_matrix[i] == NULL)
500  {
501  VecDuplicate(eq_data_->ls[i]->get_solution(), &mass_vec[i]);
502  MatMult(*(eq_data_->ls_dt[i]->get_matrix()), eq_data_->ls[i]->get_solution(), mass_vec[i]);
503  MatConvert(*( eq_data_->ls_dt[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &mass_matrix[i]);
504  }
505  else
506  MatCopy(*( eq_data_->ls_dt[i]->get_matrix() ), mass_matrix[i], DIFFERENT_NONZERO_PATTERN);
507  }
508  }
509 
510  // assemble stiffness matrix
511  if (stiffness_matrix[0] == NULL
512  || eq_fields_->subset(FieldFlag::in_main_matrix).changed()
513  || eq_fields_->flow_flux.changed())
514  {
515  // new fluxes can change the location of Neumann boundary,
516  // thus stiffness matrix must be reassembled
517  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
518  {
519  eq_data_->ls[i]->start_add_assembly();
520  eq_data_->ls[i]->mat_zero_entries();
521  }
522  stiffness_assembly_->assemble(eq_data_->dh_);
523  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
524  {
525  eq_data_->ls[i]->finish_assembly();
526 
527  if (stiffness_matrix[i] == NULL)
528  MatConvert(*( eq_data_->ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
529  else
530  MatCopy(*( eq_data_->ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
531  }
532  }
533 
534  // assemble right hand side (due to sources and boundary conditions)
535  if (rhs[0] == NULL
536  || eq_fields_->subset(FieldFlag::in_rhs).changed()
537  || eq_fields_->flow_flux.changed())
538  {
539  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
540  {
541  eq_data_->ls[i]->start_add_assembly();
542  eq_data_->ls[i]->rhs_zero_entries();
543  }
544  sources_assembly_->assemble(eq_data_->dh_);
545  bdr_cond_assembly_->assemble(eq_data_->dh_);
546  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
547  {
548  eq_data_->ls[i]->finish_assembly();
549 
550  if (rhs[i] == nullptr) VecDuplicate(*( eq_data_->ls[i]->get_rhs() ), &rhs[i]);
551  VecCopy(*( eq_data_->ls[i]->get_rhs() ), rhs[i]);
552  }
553  }
554 
555  /* Apply backward Euler time integration.
556  *
557  * Denoting A the stiffness matrix and M the mass matrix, the algebraic system at the k-th time level reads
558  *
559  * (1/dt M + A)u^k = f + 1/dt M.u^{k-1}
560  *
561  * Hence we modify at each time level the right hand side:
562  *
563  * f^k = f + 1/dt M u^{k-1},
564  *
565  * where f stands for the term stemming from the force and boundary conditions.
566  * Accordingly, we set
567  *
568  * A^k = A + 1/dt M.
569  *
570  */
571  Mat m;
572  START_TIMER("solve");
573  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
574  {
575  MatConvert(stiffness_matrix[i], MATSAME, MAT_INITIAL_MATRIX, &m);
576  MatAXPY(m, 1./Model::time_->dt(), mass_matrix[i], SUBSET_NONZERO_PATTERN);
577  eq_data_->ls[i]->set_matrix(m, DIFFERENT_NONZERO_PATTERN);
578  Vec w;
579  VecDuplicate(rhs[i], &w);
580  VecWAXPY(w, 1./Model::time_->dt(), mass_vec[i], rhs[i]);
581  eq_data_->ls[i]->set_rhs(w);
582 
583  VecDestroy(&w);
584  chkerr(MatDestroy(&m));
585 
586  eq_data_->ls[i]->solve();
587 
588  // update mass_vec due to possible changes in mass matrix
589  MatMult(*(eq_data_->ls_dt[i]->get_matrix()), eq_data_->ls[i]->get_solution(), mass_vec[i]);
590  }
591  END_TIMER("solve");
592 
593  // Possibly output matrices for debug reasons.
594  // for (unsigned int i=0; i<eq_data_->n_substances(); i++){
595  // string conc_name = eq_data_->substances().names()[i] + "_" + std::to_string(eq_data_->time_->step().index());
596  // eq_data_->ls[i]->view("stiff_" + conc_name);
597  // eq_data_->ls_dt[i]->view("mass_" + conc_name);
598  // }
599 
600  calculate_cumulative_balance();
601 
602  END_TIMER("DG-ONE STEP");
603 }
604 
605 
606 template<class Model>
608 {
609  // calculate element averages of solution
610  for (auto cell : eq_data_->dh_->own_range() )
611  {
612  LocDofVec loc_dof_indices = cell.get_loc_dof_indices();
613  unsigned int n_dofs=loc_dof_indices.n_rows;
614 
615  DHCellAccessor dh_p0_cell = eq_data_->dh_p0->cell_accessor_from_element(cell.elm_idx());
616  IntIdx dof_p0 = dh_p0_cell.get_loc_dof_indices()[0];
617 
618  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
619  {
620  eq_data_->conc_fe[sbi]->vec().set(dof_p0, 0);
621 
622  for (unsigned int j=0; j<n_dofs; ++j)
623  eq_data_->conc_fe[sbi]->vec().add( dof_p0, eq_data_->ls[sbi]->get_solution_array()[loc_dof_indices[j]] );
624 
625  eq_data_->conc_fe[sbi]->vec().set( dof_p0, max(eq_data_->conc_fe[sbi]->vec().get(dof_p0)/n_dofs, 0.) );
626  }
627  }
628 }
629 
630 
631 
632 
633 template<class Model>
635 {
636  //if (!Model::time_->is_current( Model::time_->marks().type_output() )) return;
637 
638 
639  START_TIMER("DG-OUTPUT");
640 
641  // gather the solution from all processors
642  eq_fields_->output_fields.set_time( this->time().step(), LimitSide::left);
643  //if (eq_fields_->output_fields.is_field_output_time(eq_fields_->output_field, this->time().step()) )
644  eq_fields_->output_fields.output(this->time().step());
645 
646  Model::output_data();
647 
648  START_TIMER("TOS-balance");
649  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
650  Model::balance_->calculate_instant(eq_data_->subst_idx_[sbi], eq_data_->ls[sbi]->get_solution());
651  Model::balance_->output();
652  END_TIMER("TOS-balance");
653 
654  END_TIMER("DG-OUTPUT");
655 }
656 
657 
658 template<class Model>
660 {
661  if (Model::balance_->cumulative())
662  {
663  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
664  {
665  Model::balance_->calculate_cumulative(eq_data_->subst_idx_[sbi], eq_data_->ls[sbi]->get_solution());
666 
667  // update source increment due to retardation
668  VecDot(eq_data_->ret_vec[sbi], eq_data_->ls[sbi]->get_solution(), &ret_sources[sbi]);
669 
670  Model::balance_->add_cumulative_source(eq_data_->subst_idx_[sbi], (ret_sources[sbi]-ret_sources_prev[sbi])/Model::time_->dt());
671  ret_sources_prev[sbi] = ret_sources[sbi];
672  }
673  }
674 }
675 
676 
677 
678 
679 template<class Model>
681 {
682  if(init_projection)
683  {
684  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
685  eq_data_->ls[sbi]->start_allocation();
686 
687  init_assembly_->assemble(eq_data_->dh_);
688 
689  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
690  eq_data_->ls[sbi]->start_add_assembly();
691 
692  init_assembly_->assemble(eq_data_->dh_);
693 
694  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
695  {
696  eq_data_->ls[sbi]->finish_assembly();
697  eq_data_->ls[sbi]->solve();
698  }
699  }
700  else
701  init_assembly_->assemble(eq_data_->dh_);
702 }
703 
704 
705 template<class Model>
707 {
708  if (solution_changed)
709  {
710  for (auto cell : eq_data_->dh_->own_range() )
711  {
712  LocDofVec loc_dof_indices = cell.get_loc_dof_indices();
713  unsigned int n_dofs=loc_dof_indices.n_rows;
714 
715  DHCellAccessor dh_p0_cell = eq_data_->dh_p0->cell_accessor_from_element(cell.elm_idx());
716  IntIdx dof_p0 = dh_p0_cell.get_loc_dof_indices()[0];
717 
718  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
719  {
720  double old_average = 0;
721  for (unsigned int j=0; j<n_dofs; ++j)
722  old_average += eq_data_->ls[sbi]->get_solution_array()[loc_dof_indices[j]];
723  old_average /= n_dofs;
724 
725  for (unsigned int j=0; j<n_dofs; ++j)
726  eq_data_->ls[sbi]->get_solution_array()[loc_dof_indices[j]]
727  += eq_data_->conc_fe[sbi]->vec().get(dof_p0) - old_average;
728  }
729  }
730  }
731  // update mass_vec for the case that mass matrix changes in next time step
732  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
733  MatMult(*(eq_data_->ls_dt[sbi]->get_matrix()), eq_data_->ls[sbi]->get_solution(), mass_vec[sbi]);
734 }
735 
736 
737 
738 
739 
740 
742 template class TransportDG<ConcDispersionModel>;
743 template class TransportDG<HeatTransferModel>;
744 
745 
746 
747 
FieldCommon::units
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
Definition: field_common.hh:153
EquationBase::user_fields_template
static Input::Type::Record & user_fields_template(std::string equation_name)
Template Record with common key user_fields for derived equations.
Definition: equation.cc:46
FieldFlag::in_rhs
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
LocDofVec
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
Input::Type::Bool
Class for declaration of the input of type Bool.
Definition: type_base.hh:452
fe_rt.hh
Definitions of Raviart-Thomas finite elements.
UnitSI::dimensionless
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
TransportDG::EqFields::dg_penalty
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
Definition: transport_dg.hh:153
factory.hh
LinSys
Definition: la_linsys_new.hh:169
transport_dg.hh
Discontinuous Galerkin method for equation of transport with dispersion.
EquationBase::eq_fieldset_
std::shared_ptr< FieldSet > eq_fieldset_
Definition: equation.hh:249
TransportDG::input_rec
Input::Record input_rec
Array for storing the output solution data.
Definition: transport_dg.hh:392
Input::Type::Integer
Class for declaration of the integral input data.
Definition: type_base.hh:483
Input::Record::val
const Ret val(const string &key) const
Definition: accessors_impl.hh:31
IntIdx
int IntIdx
Definition: index_types.hh:25
EquationBase::time
TimeGovernor & time()
Definition: equation.hh:151
Input::Type::Selection::close
const Selection & close() const
Close the Selection, no more values can be added.
Definition: type_selection.cc:65
fe_values.hh
Class FEValues calculates finite element data on the actual cells such as shape function values,...
value
static constexpr bool value
Definition: json.hpp:87
chkerr
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:142
FLOW123D_FORCE_LINK_IN_CHILD
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:104
LinSys_PETSC
Definition: linsys_PETSC.hh:43
TransportDG::EqFields::EqFields
EqFields()
Definition: transport_dg.cc:113
Side::dim
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
Definition: accessors_impl.hh:198
std::vector
Definition: doxy_dummy_defs.hh:7
TransportDG::update_solution
void update_solution() override
Computes the solution in one time instant.
Definition: transport_dg.cc:472
ElementAccessor< 3 >
arma::vec3
Definition: doxy_dummy_defs.hh:17
FieldCommon::flags
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Definition: field_common.hh:192
advection_diffusion_model.hh
Discontinuous Galerkin method for equation of transport with dispersion.
field_fe.hh
linsys_PETSC.hh
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
EquationBase::init_user_fields
void init_user_fields(Input::Array user_fields, FieldSet &output_fields)
Definition: equation.cc:84
TransportDG::set_initial_condition
void set_initial_condition()
Calculates the dispersivity (diffusivity) tensor from the velocity field.
Definition: transport_dg.cc:680
TransportDG::EqData::set_DG_parameters_boundary
void set_DG_parameters_boundary(Side side, const int K_size, const std::vector< arma::mat33 > &K, const double flux, const arma::vec3 &normal_vector, const double alpha, double &gamma)
Sets up parameters of the DG method on a given boundary edge.
Definition: transport_dg.cc:170
index_types.hh
TransportDG::mass_assembly_
GenericAssembly< MassAssemblyDim > * mass_assembly_
general assembly objects, hold assembly objects of appropriate dimension
Definition: transport_dg.hh:421
TransportDG::calculate_cumulative_balance
void calculate_cumulative_balance()
Definition: transport_dg.cc:659
TransportDG::compute_p0_interpolation
void compute_p0_interpolation()
Compute P0 interpolation of the solution (used in reaction term).
Definition: transport_dg.cc:607
FieldFlag::equation_external_output
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
TransportDG::EqFields::region_id
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport_dg.hh:154
TransportDG::TransportDG
TransportDG(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
Definition: transport_dg.cc:206
fe_p.hh
Definitions of basic Lagrangean finite elements with polynomial shape functions.
FieldCommon::flags_add
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
Definition: field_common.hh:198
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
equation_output.hh
EquationBase::balance
std::shared_ptr< Balance > balance() const
Definition: equation.hh:189
TransportDG::output_data
void output_data()
Postprocesses the solution and writes to output file.
Definition: transport_dg.cc:634
LimitSide::left
@ left
dh_cell_accessor.hh
accessors.hh
TransportDG::eq_fields_
std::shared_ptr< EqFields > eq_fields_
Fields for model parameters.
Definition: transport_dg.hh:359
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
sys_profiler.hh
generic_assembly.hh
quadrature_lib.hh
Definitions of particular quadrature rules on simplices.
TransportDG::initialize
void initialize() override
Definition: transport_dg.cc:246
output_time.hh
TransportDG::get_dg_variant_selection_input_type
static const Input::Type::Selection & get_dg_variant_selection_input_type()
Input type for the DG variant selection.
Definition: transport_dg.cc:54
TransportDG::preallocate
void preallocate()
Definition: transport_dg.cc:441
Side
Definition: accessors.hh:390
GenericField::subdomain
static auto subdomain(Mesh &mesh) -> IndexField
Definition: generic_field.impl.hh:58
TransportOperatorSplitting::output_data
void output_data() override
Write computed fields.
Definition: transport_operator_splitting.cc:242
Input::Type::Default::obligatory
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
FieldCommon::print_message_table
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
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
Side::element
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
Definition: accessors_impl.hh:218
TransportDG::allocation_done
bool allocation_done
Indicates whether matrices have been preallocated.
Definition: transport_dg.hh:415
TransportDG::EqFields
Definition: transport_dg.hh:147
Input::Type::Selection
Template for classes storing finite set of named values.
Definition: type_selection.hh:65
Side::n_nodes
unsigned int n_nodes() const
Returns number of nodes of the side.
Definition: accessors.hh:436
ElementAccessor::node
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:230
Input::Type::Record::close
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
heat_model.hh
Discontinuous Galerkin method for equation of transport with dispersion.
generic_field.hh
Input::Type
Definition: balance.hh:41
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
TransportDG
Transport with dispersion implemented using discontinuous Galerkin method.
Definition: transport_dg.hh:134
TransportDG::EqFields::fracture_sigma
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
Definition: transport_dg.hh:152
GenericField::region_id
static auto region_id(Mesh &mesh) -> IndexField
Definition: generic_field.impl.hh:39
FieldFlag::in_main_matrix
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:49
FieldCommon::input_default
FieldCommon & input_default(const string &input_default)
Definition: field_common.hh:140
FieldFE
Definition: field.hh:60
conc_dispersion_model.hh
TransportDG::eq_data_
std::shared_ptr< EqData > eq_data_
Data for model parameters.
Definition: transport_dg.hh:362
TransportDG::DGVariant
DGVariant
Definition: transport_dg.hh:244
Mesh
Definition: mesh.h:362
Input::Type::Record::copy_keys
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
multi_field.hh
Input::Type::Array
Class for declaration of inputs sequences.
Definition: type_base.hh:339
TransportDG::update_after_reactions
void update_after_reactions(bool solution_changed)
Definition: transport_dg.cc:706
Model
Definition: field_model.hh:338
DHCellAccessor
Cell accessor allow iterate over DOF handler cells.
Definition: dh_cell_accessor.hh:43
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
fe_value_handler.hh
TransportDG::EqFields::subdomain
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport_dg.hh:155
OutputTime::CORNER_DATA
@ CORNER_DATA
Definition: output_time.hh:110
MixedPtr< FE_P_disc >
DHCellAccessor::get_loc_dof_indices
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
Definition: dh_cell_accessor.hh:88
assembly_dg.hh
concentration_model.hh
Discontinuous Galerkin method for equation of transport with dispersion.
TransportDG::zero_time_step
void zero_time_step() override
Initialize solution in the zero time.
Definition: transport_dg.cc:408
TransportDG::~TransportDG
~TransportDG() override
Destructor.
Definition: transport_dg.cc:356
balance.hh
FieldCommon::description
FieldCommon & description(const string &description)
Definition: field_common.hh:128
TransportDG::get_input_type
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: transport_dg.cc:80
Side::node
NodeAccessor< 3 > node(unsigned int i) const
Returns node for given local index i on the side.
Definition: accessors_impl.hh:212
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
FieldFlag::in_time_term
static constexpr Mask in_time_term
A field is part of time term of the equation.
Definition: field_flag.hh:47
START_TIMER
#define START_TIMER(tag)
Starts a timer with specified tag.
Definition: sys_profiler.hh:115
TransportDG::EqData::elem_anisotropy
double elem_anisotropy(ElementAccessor< 3 > e) const
Compute and return anisotropy of given element.
Definition: transport_dg.cc:154
GenericAssembly
Generic class of assemblation.
Definition: generic_assembly.hh:165
TransportDG::EqFields::output_fields
EquationOutput output_fields
Definition: transport_dg.hh:157
LinSys_PETSC::get_input_type
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:32
END_TIMER
#define END_TIMER(tag)
Ends a timer with specified tag.
Definition: sys_profiler.hh:149
Element::n_nodes
unsigned int n_nodes() const
Definition: elements.h:125
FieldCommon::name
FieldCommon & name(const string &name)
Definition: field_common.hh:121
FieldValue
Definition: field_values.hh:645