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