Flow123d  JS_before_hm-1601-gc6ac32d
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>
53 const Selection & TransportDG<Model>::get_dg_variant_selection_input_type() {
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")
83  "Solver for the linear system.")
84  .declare_key("input_fields", Array(
86  .make_field_descriptor_type(equation_name)),
88  "Input fields of the equation.")
90  "Variant of the interior penalty discontinuous Galerkin method.")
91  .declare_key("dg_order", Integer(0,3), Default("1"),
92  "Polynomial order for the finite element in DG method (order 0 is suitable if there is no diffusion/dispersion).")
93  .declare_key("output",
94  EqFields().output_fields.make_output_type(equation_name, ""),
95  IT::Default("{ \"fields\": [ " + Model::ModelEqData::default_output_field() + "] }"),
96  "Specification of output fields and output times.")
97  .close();
98 }
99 
100 template<class Model>
102  Input::register_class< TransportDG<Model>, Mesh &, const Input::Record>(std::string(Model::ModelEqData::name()) + "_DG") +
104 
105 
106 
107 template<class Model>
109 {
110  *this+=fracture_sigma
111  .name("fracture_sigma")
112  .description(
113  "Coefficient of diffusive transfer through fractures (for each substance).")
115  .input_default("1.0")
117 
118  *this+=dg_penalty
119  .name("dg_penalty")
120  .description(
121  "Penalty parameter influencing the discontinuity of the solution (for each substance). "
122  "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
124  .input_default("1.0")
126 
127  *this += region_id.name("region_id")
130  .description("Region ids.");
131 
132  *this += subdomain.name("subdomain")
135  .description("Subdomain ids of the domain decomposition.");
136 
137 
138  // add all input fields to the output list
139  output_fields += *this;
140 
141 }
142 
143 
144 
145 // return the ratio of longest and shortest edge
146 template<class Model>
148 {
149  double h_max = 0, h_min = numeric_limits<double>::infinity();
150  for (unsigned int i=0; i<e->n_nodes(); i++)
151  for (unsigned int j=i+1; j<e->n_nodes(); j++)
152  {
153  double dist = arma::norm(*e.node(i) - *e.node(j));
154  h_max = max(h_max, dist);
155  h_min = min(h_min, dist);
156  }
157  return h_max/h_min;
158 }
159 
160 
161 
162 template<class Model>
164  const int K_size,
165  const vector<arma::mat33> &K,
166  const double flux,
167  const arma::vec3 &normal_vector,
168  const double alpha,
169  double &gamma)
170 {
171  double delta = 0, h = 0;
172 
173  // calculate the side diameter
174  if (side.dim() == 0)
175  {
176  h = 1;
177  }
178  else
179  {
180  for (unsigned int i=0; i<side.n_nodes(); i++)
181  for (unsigned int j=i+1; j<side.n_nodes(); j++) {
182  double dist = arma::norm(*side.node(i) - *side.node(j));
183  h = max(h, dist);
184  }
185 
186  }
187 
188  // delta is set to the average value of Kn.n on the side
189  for (int k=0; k<K_size; k++)
190  delta += dot(K[k]*normal_vector,normal_vector);
191  delta /= K_size;
192 
193  gamma = 0.5*fabs(flux) + alpha/h*delta*elem_anisotropy(side.element());
194 }
195 
196 
197 
198 template<typename Model>
200  : Model(init_mesh, in_rec),
201  input_rec(in_rec),
202  allocation_done(false)
203 {
204  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
205  // due to constexpr optimization.
206  START_TIMER(Model::ModelEqData::name());
207  // Check that Model is derived from AdvectionDiffusionModel.
209 
210  eq_data_ = make_shared<EqData>();
211  eq_fields_ = make_shared<EqFields>();
212  eq_fields_->add_coords_field();
213  this->eq_fieldset_ = eq_fields_.get();
214  Model::init_balance(in_rec);
215 
216 
217  // Set up physical parameters.
218  eq_fields_->set_mesh(init_mesh);
219  eq_fields_->region_id = GenericField<3>::region_id(*Model::mesh_);
220  eq_fields_->subdomain = GenericField<3>::subdomain(*Model::mesh_);
221 
222 
223  // DG data parameters
224  eq_data_->dg_variant = in_rec.val<DGVariant>("dg_variant");
225  eq_data_->dg_order = in_rec.val<unsigned int>("dg_order");
226 
227  Model::init_from_input(in_rec);
228 
229  MixedPtr<FE_P_disc> fe(eq_data_->dg_order);
230  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>(Model::mesh_, fe);
231  eq_data_->dh_ = make_shared<DOFHandlerMultiDim>(*Model::mesh_);
232  eq_data_->dh_->distribute_dofs(ds);
233  //DebugOut().fmt("TDG: solution size {}\n", eq_data_->dh_->n_global_dofs());
234 
235 }
236 
237 
238 template<class Model>
240 {
241  eq_fields_->set_components(eq_data_->substances_.names());
242  eq_fields_->set_input_list( input_rec.val<Input::Array>("input_fields"), *(Model::time_) );
243  eq_data_->set_time_governor(Model::time_);
244  eq_data_->balance_ = this->balance();
245  eq_fields_->initialize();
246 
247  // DG stabilization parameters on boundary edges
248  eq_data_->gamma.resize(eq_data_->n_substances());
249  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
250  eq_data_->gamma[sbi].resize(Model::mesh_->boundary_.size());
251 
252  // Resize coefficient arrays
253  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)));
254  ret_sources.resize(eq_data_->n_substances());
255  ret_sources_prev.resize(eq_data_->n_substances());
256 
257  output_vec.resize(eq_data_->n_substances());
258  eq_fields_->output_field.set_components(eq_data_->substances_.names());
259  eq_fields_->output_field.set_mesh(*Model::mesh_);
260  eq_fields_->output_type(OutputTime::CORNER_DATA);
261 
262  eq_fields_->output_field.setup_components();
263  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
264  {
265  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
266  auto output_field_ptr = create_field_fe< 3, FieldValue<3>::Scalar >(eq_data_->dh_);
267  eq_fields_->output_field[sbi].set(output_field_ptr, 0);
268  output_vec[sbi] = output_field_ptr->vec();
269  }
270 
271  // set time marks for writing the output
272  eq_fields_->output_fields.initialize(Model::output_stream_, Model::mesh_, input_rec.val<Input::Record>("output"), this->time());
273 
274  // equation default PETSc solver options
275  std::string petsc_default_opts;
276  if (eq_data_->dh_->distr()->np() == 1)
277  petsc_default_opts = "-ksp_type bcgs -pc_type ilu -pc_factor_levels 2 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
278  else
279  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";
280 
281  // allocate matrix and vector structures
282  eq_data_->ls = new LinSys*[eq_data_->n_substances()];
283  eq_data_->ls_dt = new LinSys*[eq_data_->n_substances()];
284  eq_data_->conc_fe.resize(eq_data_->n_substances());
285 
286  MixedPtr<FE_P_disc> fe(0);
287  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>(Model::mesh_, fe);
288  eq_data_->dh_p0 = make_shared<DOFHandlerMultiDim>(*Model::mesh_);
289  eq_data_->dh_p0->distribute_dofs(ds);
290 
291  stiffness_matrix.resize(eq_data_->n_substances(), nullptr);
292  mass_matrix.resize(eq_data_->n_substances(), nullptr);
293  rhs.resize(eq_data_->n_substances(), nullptr);
294  mass_vec.resize(eq_data_->n_substances(), nullptr);
295  eq_data_->ret_vec.resize(eq_data_->n_substances(), nullptr);
296 
297  for (unsigned int sbi = 0; sbi < eq_data_->n_substances(); sbi++) {
298  eq_data_->ls[sbi] = new LinSys_PETSC(eq_data_->dh_->distr().get(), petsc_default_opts);
299  ( (LinSys_PETSC *)eq_data_->ls[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
300  eq_data_->ls[sbi]->set_solution(output_vec[sbi].petsc_vec());
301 
302  eq_data_->ls_dt[sbi] = new LinSys_PETSC(eq_data_->dh_->distr().get(), petsc_default_opts);
303  ( (LinSys_PETSC *)eq_data_->ls_dt[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
304 
305  eq_data_->conc_fe[sbi] = create_field_fe< 3, FieldValue<3>::Scalar >(eq_data_->dh_p0);
306 
307  VecDuplicate(eq_data_->ls[sbi]->get_solution(), &eq_data_->ret_vec[sbi]);
308  }
309 
310 
311  // create assemblation object, finite element structures and distribute DOFs
317 
318  // initialization of balance object
319  Model::balance_->allocate(eq_data_->dh_->distr()->lsize(), mass_assembly_->eval_points()->max_size());
320 
321  int qsize = mass_assembly_->eval_points()->max_size();
322  eq_data_->dif_coef.resize(eq_data_->n_substances());
323  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
324  {
325  eq_data_->dif_coef[sbi].resize(qsize);
326  }
327 }
328 
329 
330 template<class Model>
332 {
333  delete Model::time_;
334 
335  if (eq_data_->gamma.size() > 0) {
336  // initialize called
337 
338  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
339  {
340  delete eq_data_->ls[i];
341  delete eq_data_->ls_dt[i];
342 
343  if (stiffness_matrix[i])
344  chkerr(MatDestroy(&stiffness_matrix[i]));
345  if (mass_matrix[i])
346  chkerr(MatDestroy(&mass_matrix[i]));
347  if (rhs[i])
348  chkerr(VecDestroy(&rhs[i]));
349  if (mass_vec[i])
350  chkerr(VecDestroy(&mass_vec[i]));
351  if (eq_data_->ret_vec[i])
352  chkerr(VecDestroy(&eq_data_->ret_vec[i]));
353  }
354  delete[] eq_data_->ls;
355  delete[] eq_data_->ls_dt;
356  //delete[] stiffness_matrix;
357  //delete[] mass_matrix;
358  //delete[] rhs;
359  //delete[] mass_vec;
360  //delete[] ret_vec;
361 
362  delete mass_assembly_;
363  delete stiffness_assembly_;
364  delete sources_assembly_;
365  delete bdr_cond_assembly_;
366  delete init_cond_assembly_;
367  }
368 
369 }
370 
371 
372 template<class Model>
374 {
375  START_TIMER(Model::ModelEqData::name());
376  eq_fields_->mark_input_times( *(Model::time_) );
377  eq_fields_->set_time(Model::time_->step(), LimitSide::left);
378  std::stringstream ss; // print warning message with table of uninitialized fields
379  if ( FieldCommon::print_message_table(ss, "transport DG") ) {
380  WarningOut() << ss.str();
381  }
382 
383 
384  // set initial conditions
386  for (unsigned int sbi = 0; sbi < eq_data_->n_substances(); sbi++)
387  ( (LinSys_PETSC *)eq_data_->ls[sbi] )->set_initial_guess_nonzero();
388 
389  // check first time assembly - needs preallocation
391 
392  // after preallocation we assemble the matrices and vectors required for mass balance
393  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
394  {
395  Model::balance_->calculate_instant(eq_data_->subst_idx_[sbi], eq_data_->ls[sbi]->get_solution());
396 
397  // add sources due to sorption
398  ret_sources_prev[sbi] = 0;
399  }
400 
401  output_data();
402 }
403 
404 
405 template<class Model>
407 {
408  // preallocate system matrix
409  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
410  {
411  // preallocate system matrix
412  eq_data_->ls[i]->start_allocation();
413  stiffness_matrix[i] = NULL;
414  rhs[i] = NULL;
415 
416  // preallocate mass matrix
417  eq_data_->ls_dt[i]->start_allocation();
418  mass_matrix[i] = NULL;
419  VecZeroEntries(eq_data_->ret_vec[i]);
420  }
421  stiffness_assembly_->assemble(eq_data_->dh_);
422  mass_assembly_->assemble(eq_data_->dh_);
423  sources_assembly_->assemble(eq_data_->dh_);
424  bdr_cond_assembly_->assemble(eq_data_->dh_);
425  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
426  {
427  VecAssemblyBegin(eq_data_->ret_vec[i]);
428  VecAssemblyEnd(eq_data_->ret_vec[i]);
429  }
430 
431  allocation_done = true;
432 }
433 
434 
435 
436 template<class Model>
438 {
439  START_TIMER("DG-ONE STEP");
440 
441  Model::time_->next_time();
442  Model::time_->view("TDG");
443 
444  START_TIMER("data reinit");
445  eq_fields_->set_time(Model::time_->step(), LimitSide::left);
446  END_TIMER("data reinit");
447 
448  // assemble mass matrix
449  if (mass_matrix[0] == NULL || eq_fields_->subset(FieldFlag::in_time_term).changed() )
450  {
451  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
452  {
453  eq_data_->ls_dt[i]->start_add_assembly();
454  eq_data_->ls_dt[i]->mat_zero_entries();
455  VecZeroEntries(eq_data_->ret_vec[i]);
456  }
457  mass_assembly_->assemble(eq_data_->dh_);
458  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
459  {
460  eq_data_->ls_dt[i]->finish_assembly();
461  VecAssemblyBegin(eq_data_->ret_vec[i]);
462  VecAssemblyEnd(eq_data_->ret_vec[i]);
463  // construct mass_vec for initial time
464  if (mass_matrix[i] == NULL)
465  {
466  VecDuplicate(eq_data_->ls[i]->get_solution(), &mass_vec[i]);
467  MatMult(*(eq_data_->ls_dt[i]->get_matrix()), eq_data_->ls[i]->get_solution(), mass_vec[i]);
468  MatConvert(*( eq_data_->ls_dt[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &mass_matrix[i]);
469  }
470  else
471  MatCopy(*( eq_data_->ls_dt[i]->get_matrix() ), mass_matrix[i], DIFFERENT_NONZERO_PATTERN);
472  }
473  }
474 
475  // assemble stiffness matrix
476  if (stiffness_matrix[0] == NULL
477  || eq_fields_->subset(FieldFlag::in_main_matrix).changed()
478  || eq_fields_->flow_flux.changed())
479  {
480  // new fluxes can change the location of Neumann boundary,
481  // thus stiffness matrix must be reassembled
482  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
483  {
484  eq_data_->ls[i]->start_add_assembly();
485  eq_data_->ls[i]->mat_zero_entries();
486  }
487  stiffness_assembly_->assemble(eq_data_->dh_);
488  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
489  {
490  eq_data_->ls[i]->finish_assembly();
491 
492  if (stiffness_matrix[i] == NULL)
493  MatConvert(*( eq_data_->ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
494  else
495  MatCopy(*( eq_data_->ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
496  }
497  }
498 
499  // assemble right hand side (due to sources and boundary conditions)
500  if (rhs[0] == NULL
501  || eq_fields_->subset(FieldFlag::in_rhs).changed()
502  || eq_fields_->flow_flux.changed())
503  {
504  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
505  {
506  eq_data_->ls[i]->start_add_assembly();
507  eq_data_->ls[i]->rhs_zero_entries();
508  }
509  sources_assembly_->assemble(eq_data_->dh_);
510  bdr_cond_assembly_->assemble(eq_data_->dh_);
511  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
512  {
513  eq_data_->ls[i]->finish_assembly();
514 
515  if (rhs[i] == nullptr) VecDuplicate(*( eq_data_->ls[i]->get_rhs() ), &rhs[i]);
516  VecCopy(*( eq_data_->ls[i]->get_rhs() ), rhs[i]);
517  }
518  }
519 
520 
521  /* Apply backward Euler time integration.
522  *
523  * Denoting A the stiffness matrix and M the mass matrix, the algebraic system at the k-th time level reads
524  *
525  * (1/dt M + A)u^k = f + 1/dt M.u^{k-1}
526  *
527  * Hence we modify at each time level the right hand side:
528  *
529  * f^k = f + 1/dt M u^{k-1},
530  *
531  * where f stands for the term stemming from the force and boundary conditions.
532  * Accordingly, we set
533  *
534  * A^k = A + 1/dt M.
535  *
536  */
537  Mat m;
538  START_TIMER("solve");
539  for (unsigned int i=0; i<eq_data_->n_substances(); i++)
540  {
541  MatConvert(stiffness_matrix[i], MATSAME, MAT_INITIAL_MATRIX, &m);
542  MatAXPY(m, 1./Model::time_->dt(), mass_matrix[i], SUBSET_NONZERO_PATTERN);
543  eq_data_->ls[i]->set_matrix(m, DIFFERENT_NONZERO_PATTERN);
544  Vec w;
545  VecDuplicate(rhs[i], &w);
546  VecWAXPY(w, 1./Model::time_->dt(), mass_vec[i], rhs[i]);
547  eq_data_->ls[i]->set_rhs(w);
548 
549  VecDestroy(&w);
550  chkerr(MatDestroy(&m));
551 
552  eq_data_->ls[i]->solve();
553 
554  // update mass_vec due to possible changes in mass matrix
555  MatMult(*(eq_data_->ls_dt[i]->get_matrix()), eq_data_->ls[i]->get_solution(), mass_vec[i]);
556  }
557  END_TIMER("solve");
558 
560 
561  END_TIMER("DG-ONE STEP");
562 }
563 
564 
565 template<class Model>
567 {
568  // calculate element averages of solution
569  for (auto cell : eq_data_->dh_->own_range() )
570  {
571  LocDofVec loc_dof_indices = cell.get_loc_dof_indices();
572  unsigned int n_dofs=loc_dof_indices.n_rows;
573 
574  DHCellAccessor dh_p0_cell = eq_data_->dh_p0->cell_accessor_from_element(cell.elm_idx());
575  IntIdx dof_p0 = dh_p0_cell.get_loc_dof_indices()[0];
576 
577  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
578  {
579  eq_data_->conc_fe[sbi]->vec().set(dof_p0, 0);
580 
581  for (unsigned int j=0; j<n_dofs; ++j)
582  eq_data_->conc_fe[sbi]->vec().add( dof_p0, eq_data_->ls[sbi]->get_solution_array()[loc_dof_indices[j]] );
583 
584  eq_data_->conc_fe[sbi]->vec().set( dof_p0, max(eq_data_->conc_fe[sbi]->vec().get(dof_p0)/n_dofs, 0.) );
585  }
586  }
587 }
588 
589 
590 
591 
592 template<class Model>
594 {
595  //if (!Model::time_->is_current( Model::time_->marks().type_output() )) return;
596 
597 
598  START_TIMER("DG-OUTPUT");
599 
600  // gather the solution from all processors
601  eq_fields_->output_fields.set_time( this->time().step(), LimitSide::left);
602  //if (eq_fields_->output_fields.is_field_output_time(eq_fields_->output_field, this->time().step()) )
603  eq_fields_->output_fields.output(this->time().step());
604 
605  Model::output_data();
606 
607  START_TIMER("TOS-balance");
608  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
609  Model::balance_->calculate_instant(eq_data_->subst_idx_[sbi], eq_data_->ls[sbi]->get_solution());
610  Model::balance_->output();
611  END_TIMER("TOS-balance");
612 
613  END_TIMER("DG-OUTPUT");
614 }
615 
616 
617 template<class Model>
619 {
620  if (Model::balance_->cumulative())
621  {
622  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
623  {
624  Model::balance_->calculate_cumulative(eq_data_->subst_idx_[sbi], eq_data_->ls[sbi]->get_solution());
625 
626  // update source increment due to retardation
627  VecDot(eq_data_->ret_vec[sbi], eq_data_->ls[sbi]->get_solution(), &ret_sources[sbi]);
628 
629  Model::balance_->add_cumulative_source(eq_data_->subst_idx_[sbi], (ret_sources[sbi]-ret_sources_prev[sbi])/Model::time_->dt());
630  ret_sources_prev[sbi] = ret_sources[sbi];
631  }
632  }
633 }
634 
635 
636 
637 
638 template<class Model>
640 {
641  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
642  eq_data_->ls[sbi]->start_allocation();
643  init_cond_assembly_->assemble(eq_data_->dh_);
644 
645  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
646  eq_data_->ls[sbi]->start_add_assembly();
647  init_cond_assembly_->assemble(eq_data_->dh_);
648 
649  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
650  {
651  eq_data_->ls[sbi]->finish_assembly();
652  eq_data_->ls[sbi]->solve();
653  }
654 }
655 
656 
657 template<class Model>
659 {
660  el_4_loc = Model::mesh_->get_el_4_loc();
661  el_ds = Model::mesh_->get_el_ds();
662 }
663 
664 
665 template<class Model>
667 {
668  if (solution_changed)
669  {
670  for (auto cell : eq_data_->dh_->own_range() )
671  {
672  LocDofVec loc_dof_indices = cell.get_loc_dof_indices();
673  unsigned int n_dofs=loc_dof_indices.n_rows;
674 
675  DHCellAccessor dh_p0_cell = eq_data_->dh_p0->cell_accessor_from_element(cell.elm_idx());
676  IntIdx dof_p0 = dh_p0_cell.get_loc_dof_indices()[0];
677 
678  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
679  {
680  double old_average = 0;
681  for (unsigned int j=0; j<n_dofs; ++j)
682  old_average += eq_data_->ls[sbi]->get_solution_array()[loc_dof_indices[j]];
683  old_average /= n_dofs;
684 
685  for (unsigned int j=0; j<n_dofs; ++j)
686  eq_data_->ls[sbi]->get_solution_array()[loc_dof_indices[j]]
687  += eq_data_->conc_fe[sbi]->vec().get(dof_p0) - old_average;
688  }
689  }
690  }
691  // update mass_vec for the case that mass matrix changes in next time step
692  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
693  MatMult(*(eq_data_->ls_dt[sbi]->get_matrix()), eq_data_->ls[sbi]->get_solution(), mass_vec[sbi]);
694 }
695 
696 template<class Model>
698 {
699  return Model::mesh_->get_row_4_el();
700 }
701 
702 
703 
704 
705 
706 
708 template class TransportDG<HeatTransferModel>;
709 
710 
711 
712 
TimeGovernor & time()
Definition: equation.hh:149
Input::Record input_rec
Record with input specification.
GenericAssembly< SourcesAssemblyDim > * sources_assembly_
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
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:566
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
std::shared_ptr< EvalPoints > eval_points() const
Geter to EvalPoints object.
double elem_anisotropy(ElementAccessor< 3 > e) const
Compute and return anisotropy of given element.
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...
Transport with dispersion implemented using discontinuous Galerkin method.
GenericAssembly< StiffnessAssemblyDim > * stiffness_assembly_
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
int IntIdx
Definition: index_types.hh:25
GenericAssembly< MassAssemblyDim > * mass_assembly_
general assembly objects, hold assembly objects of appropriate dimension
Abstract linear system class.
Definition: balance.hh:40
vector< VectorMPI > output_vec
Array for storing the output solution data.
Field< 3, FieldValue< 3 >::Scalar > region_id
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:77
Fields computed from the mesh data.
void set_initial_condition()
Calculates the dispersivity (diffusivity) tensor from the velocity field.
Cell accessor allow iterate over DOF handler cells.
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
std::shared_ptr< EqFields > eq_fields_
Fields for model parameters.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
Discontinuous Galerkin method for equation of transport with dispersion.
Class for declaration of the integral input data.
Definition: type_base.hh:483
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
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:200
Class for declaration of inputs sequences.
Definition: type_base.hh:339
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.
void preallocate()
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
void compute_p0_interpolation()
Compute P0 interpolation of the solution (used in reaction term).
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
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
const Ret val(const string &key) const
GenericAssembly< BdrConditionAssemblyDim > * bdr_cond_assembly_
static auto region_id(Mesh &mesh) -> IndexField
FieldSet * eq_fieldset_
Definition: equation.hh:227
#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
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
Field< 3, FieldValue< 3 >::Scalar > subdomain
GenericAssembly< InitConditionAssemblyDim > * init_cond_assembly_
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.
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.
std::shared_ptr< EqData > eq_data_
Data for model parameters.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: transport_dg.cc:79
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.
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
std::shared_ptr< Balance > balance() const
Access to balance object of Model.
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:270
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.
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.
unsigned int n_nodes() const
Returns number of nodes of the side.
Definition: accessors.hh:407
EquationOutput output_fields