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