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