Flow123d  release_3.0.0-684-g928e266
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/sys_profiler.hh"
21 
22 #include "io/output_time.hh"
24 #include "fem/mapping_p1.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 "flow/mh_dofhandler.hh"
35 #include "transport/heat_model.hh"
36 #include "coupling/balance.hh"
37 
38 #include "fields/multi_field.hh"
39 #include "fields/generic_field.hh"
40 #include "input/factory.hh"
42 #include "mesh/long_idx.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")
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  EqData().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 
108 FEObjects::FEObjects(Mesh *mesh_, unsigned int fe_order)
109 {
110  unsigned int q_order;
111 
112  q_order = 2*fe_order;
113  fe0_ = new FE_P_disc<0>(fe_order);
114  fe1_ = new FE_P_disc<1>(fe_order);
115  fe2_ = new FE_P_disc<2>(fe_order);
116  fe3_ = new FE_P_disc<3>(fe_order);
117 
118  fe_rt1_ = new FE_RT0<1>;
119  fe_rt2_ = new FE_RT0<2>;
120  fe_rt3_ = new FE_RT0<3>;
121 
122  q0_ = new QGauss<0>(q_order);
123  q1_ = new QGauss<1>(q_order);
124  q2_ = new QGauss<2>(q_order);
125  q3_ = new QGauss<3>(q_order);
126 
127  map1_ = new MappingP1<1,3>;
128  map2_ = new MappingP1<2,3>;
129  map3_ = new MappingP1<3,3>;
130 
131  ds_ = std::make_shared<EqualOrderDiscreteSpace>(mesh_, fe0_, fe1_, fe2_, fe3_);
132  dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
133 
134  dh_->distribute_dofs(ds_);
135 }
136 
137 
139 {
140  delete fe0_;
141  delete fe1_;
142  delete fe2_;
143  delete fe3_;
144  delete fe_rt1_;
145  delete fe_rt2_;
146  delete fe_rt3_;
147  delete q0_;
148  delete q1_;
149  delete q2_;
150  delete q3_;
151  delete map1_;
152  delete map2_;
153  delete map3_;
154 }
155 
156 template<> FiniteElement<0> *FEObjects::fe<0>() { return fe0_; }
157 template<> FiniteElement<1> *FEObjects::fe<1>() { return fe1_; }
158 template<> FiniteElement<2> *FEObjects::fe<2>() { return fe2_; }
159 template<> FiniteElement<3> *FEObjects::fe<3>() { return fe3_; }
160 
161 template<> FiniteElement<0> *FEObjects::fe_rt<0>() { return 0; }
162 template<> FiniteElement<1> *FEObjects::fe_rt<1>() { return fe_rt1_; }
163 template<> FiniteElement<2> *FEObjects::fe_rt<2>() { return fe_rt2_; }
164 template<> FiniteElement<3> *FEObjects::fe_rt<3>() { return fe_rt3_; }
165 
166 template<> Quadrature<0> *FEObjects::q<0>() { return q0_; }
167 template<> Quadrature<1> *FEObjects::q<1>() { return q1_; }
168 template<> Quadrature<2> *FEObjects::q<2>() { return q2_; }
169 template<> Quadrature<3> *FEObjects::q<3>() { return q3_; }
170 
171 template<> MappingP1<1,3> *FEObjects::mapping<1>() { return map1_; }
172 template<> MappingP1<2,3> *FEObjects::mapping<2>() { return map2_; }
173 template<> MappingP1<3,3> *FEObjects::mapping<3>() { return map3_; }
174 
175 std::shared_ptr<DOFHandlerMultiDim> FEObjects::dh() { return dh_; }
176 
177 
178 template<class Model>
179 TransportDG<Model>::EqData::EqData() : Model::ModelEqData()
180 {
181  *this+=fracture_sigma
182  .name("fracture_sigma")
183  .description(
184  "Coefficient of diffusive transfer through fractures (for each substance).")
186  .input_default("1.0")
188 
189  *this+=dg_penalty
190  .name("dg_penalty")
191  .description(
192  "Penalty parameter influencing the discontinuity of the solution (for each substance). "
193  "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
195  .input_default("1.0")
197 
198  *this += region_id.name("region_id")
201  .description("Region ids.");
202 
203  *this += subdomain.name("subdomain")
206  .description("Subdomain ids of the domain decomposition.");
207 
208 
209  // add all input fields to the output list
210  output_fields += *this;
211 
212 }
213 
214 template<class Model>
216  : Model(init_mesh, in_rec),
217  input_rec(in_rec),
218  allocation_done(false)
219 {
220  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
221  // due to constexpr optimization.
222  START_TIMER(Model::ModelEqData::name());
223  // Check that Model is derived from AdvectionDiffusionModel.
225 
226  this->eq_data_ = &data_;
227 
228 
229  // Set up physical parameters.
230  data_.set_mesh(init_mesh);
231  data_.region_id = GenericField<3>::region_id(*Model::mesh_);
232  data_.subdomain = GenericField<3>::subdomain(*Model::mesh_);
233 
234 
235  // DG variant and order
236  dg_variant = in_rec.val<DGVariant>("dg_variant");
237  dg_order = in_rec.val<unsigned int>("dg_order");
238 
239  Model::init_from_input(in_rec);
240 
241  // create finite element structures and distribute DOFs
242  feo = new FEObjects(Model::mesh_, dg_order);
243  //DebugOut().fmt("TDG: solution size {}\n", feo->dh()->n_global_dofs());
244 
245 }
246 
247 
248 template<class Model>
250 {
251  data_.set_components(Model::substances_.names());
252  data_.set_input_list( input_rec.val<Input::Array>("input_fields"), *(Model::time_) );
253 
254  // DG stabilization parameters on boundary edges
255  gamma.resize(Model::n_substances());
256  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
257  gamma[sbi].resize(Model::mesh_->boundary_.size());
258 
259  // Resize coefficient arrays
260  int qsize = max(feo->q<0>()->size(), max(feo->q<1>()->size(), max(feo->q<2>()->size(), feo->q<3>()->size())));
261  int max_edg_sides = max(Model::mesh_->max_edge_sides(1), max(Model::mesh_->max_edge_sides(2), Model::mesh_->max_edge_sides(3)));
262  mm_coef.resize(qsize);
263  ret_coef.resize(Model::n_substances());
264  ret_sources.resize(Model::n_substances());
265  ret_sources_prev.resize(Model::n_substances());
266  ad_coef.resize(Model::n_substances());
267  dif_coef.resize(Model::n_substances());
268  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
269  {
270  ret_coef[sbi].resize(qsize);
271  ad_coef[sbi].resize(qsize);
272  dif_coef[sbi].resize(qsize);
273  }
274  ad_coef_edg.resize(max_edg_sides);
275  dif_coef_edg.resize(max_edg_sides);
276  for (int sd=0; sd<max_edg_sides; sd++)
277  {
278  ad_coef_edg[sd].resize(Model::n_substances());
279  dif_coef_edg[sd].resize(Model::n_substances());
280  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
281  {
282  ad_coef_edg[sd][sbi].resize(qsize);
283  dif_coef_edg[sd][sbi].resize(qsize);
284  }
285  }
286 
287  output_vec.resize(Model::n_substances());
288  //output_solution.resize(Model::n_substances());
289  int rank;
290  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
291  unsigned int output_vector_size= (rank==0)?feo->dh()->n_global_dofs():0;
292  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
293  {
294  // for each substance we allocate output array and vector
295  //output_solution[sbi] = new double[feo->dh()->n_global_dofs()];
296  output_vec[sbi].resize(output_vector_size);
297  }
298  data_.output_field.set_components(Model::substances_.names());
299  data_.output_field.set_mesh(*Model::mesh_);
300  data_.output_type(OutputTime::CORNER_DATA);
301 
302  data_.output_field.setup_components();
303  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
304  {
305  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
306  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar> > output_field_ptr(new FieldFE<3, FieldValue<3>::Scalar>);
307  output_field_ptr->set_fe_data(feo->dh()->sequential(), feo->mapping<1>(), feo->mapping<2>(), feo->mapping<3>(), &output_vec[sbi]);
308  data_.output_field[sbi].set_field(Model::mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
309  }
310 
311  // set time marks for writing the output
312  data_.output_fields.initialize(Model::output_stream_, Model::mesh_, input_rec.val<Input::Record>("output"), this->time());
313 
314  // equation default PETSc solver options
315  std::string petsc_default_opts;
316  if (feo->dh()->distr()->np() == 1)
317  petsc_default_opts = "-ksp_type bcgs -pc_type ilu -pc_factor_levels 2 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
318  else
319  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";
320 
321  // allocate matrix and vector structures
322  ls = new LinSys*[Model::n_substances()];
323  ls_dt = new LinSys*[Model::n_substances()];
324  solution_elem_ = new double*[Model::n_substances()];
325 
326  stiffness_matrix.resize(Model::n_substances(), nullptr);
327  mass_matrix.resize(Model::n_substances(), nullptr);
328  rhs.resize(Model::n_substances(), nullptr);
329  mass_vec.resize(Model::n_substances(), nullptr);
330  ret_vec.resize(Model::n_substances(), nullptr);
331 
332  for (unsigned int sbi = 0; sbi < Model::n_substances(); sbi++) {
333  ls[sbi] = new LinSys_PETSC(feo->dh()->distr().get(), petsc_default_opts);
334  ( (LinSys_PETSC *)ls[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
335  ls[sbi]->set_solution(NULL);
336 
337  ls_dt[sbi] = new LinSys_PETSC(feo->dh()->distr().get(), petsc_default_opts);
338  ( (LinSys_PETSC *)ls_dt[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
339  solution_elem_[sbi] = new double[Model::mesh_->get_el_ds()->lsize()];
340 
341  VecDuplicate(ls[sbi]->get_solution(), &ret_vec[sbi]);
342  }
343 
344 
345  // initialization of balance object
346  Model::balance_->allocate(feo->dh()->distr()->lsize(),
347  max(feo->fe<1>()->n_dofs(), max(feo->fe<2>()->n_dofs(), feo->fe<3>()->n_dofs())));
348 
349 }
350 
351 
352 template<class Model>
354 {
355  delete Model::time_;
356 
357  if (gamma.size() > 0) {
358  // initialize called
359 
360  for (unsigned int i=0; i<Model::n_substances(); i++)
361  {
362  delete ls[i];
363  delete[] solution_elem_[i];
364  delete ls_dt[i];
365 
366  if (stiffness_matrix[i])
367  chkerr(MatDestroy(&stiffness_matrix[i]));
368  if (mass_matrix[i])
369  chkerr(MatDestroy(&mass_matrix[i]));
370  if (rhs[i])
371  chkerr(VecDestroy(&rhs[i]));
372  if (mass_vec[i])
373  chkerr(VecDestroy(&mass_vec[i]));
374  if (ret_vec[i])
375  chkerr(VecDestroy(&ret_vec[i]));
376  }
377  delete[] ls;
378  delete[] solution_elem_;
379  delete[] ls_dt;
380  //delete[] stiffness_matrix;
381  //delete[] mass_matrix;
382  //delete[] rhs;
383  //delete[] mass_vec;
384  //delete[] ret_vec;
385  delete feo;
386 
387  }
388 
389 }
390 
391 
392 template<class Model>
394 {
395  VecScatter output_scatter;
396  VecScatterCreateToZero(ls[0]->get_solution(), &output_scatter, PETSC_NULL);
397  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
398  {
399  // gather solution to output_vec[sbi]
400  VecScatterBegin(output_scatter, ls[sbi]->get_solution(), (output_vec[sbi]).petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
401  VecScatterEnd(output_scatter, ls[sbi]->get_solution(), (output_vec[sbi]).petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
402  }
403  chkerr(VecScatterDestroy(&(output_scatter)));
404 
405 }
406 
407 
408 
409 template<class Model>
411 {
412  START_TIMER(Model::ModelEqData::name());
413  data_.mark_input_times( *(Model::time_) );
414  data_.set_time(Model::time_->step(), LimitSide::left);
415  std::stringstream ss; // print warning message with table of uninitialized fields
416  if ( FieldCommon::print_message_table(ss, "transport DG") ) {
417  WarningOut() << ss.str();
418  }
419 
420 
421  // set initial conditions
423  for (unsigned int sbi = 0; sbi < Model::n_substances(); sbi++)
424  ( (LinSys_PETSC *)ls[sbi] )->set_initial_guess_nonzero();
425 
426  // check first time assembly - needs preallocation
428 
429  // after preallocation we assemble the matrices and vectors required for mass balance
430  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
431  {
432  Model::balance_->calculate_instant(Model::subst_idx[sbi], ls[sbi]->get_solution());
433 
434  // add sources due to sorption
435  ret_sources_prev[sbi] = 0;
436  }
437 
438  output_data();
439 }
440 
441 
442 template<class Model>
444 {
445  // preallocate system matrix
446  for (unsigned int i=0; i<Model::n_substances(); i++)
447  {
448  // preallocate system matrix
449  ls[i]->start_allocation();
450  stiffness_matrix[i] = NULL;
451  rhs[i] = NULL;
452 
453  // preallocate mass matrix
454  ls_dt[i]->start_allocation();
455  mass_matrix[i] = NULL;
456  VecZeroEntries(ret_vec[i]);
457  }
460  set_sources();
462  for (unsigned int i=0; i<Model::n_substances(); i++)
463  {
464  VecAssemblyBegin(ret_vec[i]);
465  VecAssemblyEnd(ret_vec[i]);
466  }
467 
468  allocation_done = true;
469 }
470 
471 
472 
473 template<class Model>
475 {
476  START_TIMER("DG-ONE STEP");
477 
478  Model::time_->next_time();
479  Model::time_->view("TDG");
480 
481  START_TIMER("data reinit");
482  data_.set_time(Model::time_->step(), LimitSide::left);
483  END_TIMER("data reinit");
484 
485  // assemble mass matrix
486  if (mass_matrix[0] == NULL || data_.subset(FieldFlag::in_time_term).changed() )
487  {
488  for (unsigned int i=0; i<Model::n_substances(); i++)
489  {
491  ls_dt[i]->mat_zero_entries();
492  VecZeroEntries(ret_vec[i]);
493  }
495  for (unsigned int i=0; i<Model::n_substances(); i++)
496  {
497  ls_dt[i]->finish_assembly();
498  VecAssemblyBegin(ret_vec[i]);
499  VecAssemblyEnd(ret_vec[i]);
500  // construct mass_vec for initial time
501  if (mass_matrix[i] == NULL)
502  {
503  VecDuplicate(ls[i]->get_solution(), &mass_vec[i]);
504  MatMult(*(ls_dt[i]->get_matrix()), ls[i]->get_solution(), mass_vec[i]);
505  MatConvert(*( ls_dt[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &mass_matrix[i]);
506  }
507  else
508  MatCopy(*( ls_dt[i]->get_matrix() ), mass_matrix[i], DIFFERENT_NONZERO_PATTERN);
509  }
510  }
511 
512  // assemble stiffness matrix
513  if (stiffness_matrix[0] == NULL
514  || data_.subset(FieldFlag::in_main_matrix).changed()
515  || Model::flux_changed)
516  {
517  // new fluxes can change the location of Neumann boundary,
518  // thus stiffness matrix must be reassembled
519  for (unsigned int i=0; i<Model::n_substances(); i++)
520  {
521  ls[i]->start_add_assembly();
522  ls[i]->mat_zero_entries();
523  }
525  for (unsigned int i=0; i<Model::n_substances(); i++)
526  {
527  ls[i]->finish_assembly();
528 
529  if (stiffness_matrix[i] == NULL)
530  MatConvert(*( ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
531  else
532  MatCopy(*( ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
533  }
534  }
535 
536  // assemble right hand side (due to sources and boundary conditions)
537  if (rhs[0] == NULL
538  || data_.subset(FieldFlag::in_rhs).changed()
539  || Model::flux_changed)
540  {
541  for (unsigned int i=0; i<Model::n_substances(); i++)
542  {
543  ls[i]->start_add_assembly();
544  ls[i]->rhs_zero_entries();
545  }
546  set_sources();
548  for (unsigned int i=0; i<Model::n_substances(); i++)
549  {
550  ls[i]->finish_assembly();
551 
552  if (rhs[i] == nullptr) VecDuplicate(*( ls[i]->get_rhs() ), &rhs[i]);
553  VecCopy(*( ls[i]->get_rhs() ), rhs[i]);
554  }
555  }
556 
557  Model::flux_changed = false;
558 
559 
560  /* Apply backward Euler time integration.
561  *
562  * Denoting A the stiffness matrix and M the mass matrix, the algebraic system at the k-th time level reads
563  *
564  * (1/dt M + A)u^k = f + 1/dt M.u^{k-1}
565  *
566  * Hence we modify at each time level the right hand side:
567  *
568  * f^k = f + 1/dt M u^{k-1},
569  *
570  * where f stands for the term stemming from the force and boundary conditions.
571  * Accordingly, we set
572  *
573  * A^k = A + 1/dt M.
574  *
575  */
576  Mat m;
577  START_TIMER("solve");
578  for (unsigned int i=0; i<Model::n_substances(); i++)
579  {
580  MatConvert(stiffness_matrix[i], MATSAME, MAT_INITIAL_MATRIX, &m);
581  MatAXPY(m, 1./Model::time_->dt(), mass_matrix[i], SUBSET_NONZERO_PATTERN);
582  ls[i]->set_matrix(m, DIFFERENT_NONZERO_PATTERN);
583  Vec w;
584  VecDuplicate(rhs[i], &w);
585  VecWAXPY(w, 1./Model::time_->dt(), mass_vec[i], rhs[i]);
586  ls[i]->set_rhs(w);
587 
588  VecDestroy(&w);
589  chkerr(MatDestroy(&m));
590 
591  ls[i]->solve();
592 
593  // update mass_vec due to possible changes in mass matrix
594  MatMult(*(ls_dt[i]->get_matrix()), ls[i]->get_solution(), mass_vec[i]);
595  }
596  END_TIMER("solve");
597 
599 
600  END_TIMER("DG-ONE STEP");
601 }
602 
603 
604 template<class Model>
606 {
607  // calculate element averages of solution
608  unsigned int i_cell=0;
609  for (auto cell : feo->dh()->own_range() )
610  {
611 
612  unsigned int n_dofs;
613  switch (cell.dim())
614  {
615  case 1:
616  n_dofs = feo->fe<1>()->n_dofs();
617  break;
618  case 2:
619  n_dofs = feo->fe<2>()->n_dofs();
620  break;
621  case 3:
622  n_dofs = feo->fe<3>()->n_dofs();
623  break;
624  }
625 
626  std::vector<LongIdx> dof_indices(n_dofs);
627  cell.get_dof_indices(dof_indices);
628 
629  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
630  {
631  solution_elem_[sbi][i_cell] = 0;
632 
633  for (unsigned int j=0; j<n_dofs; ++j)
634  solution_elem_[sbi][i_cell] += ls[sbi]->get_solution_array()[dof_indices[j]-feo->dh()->distr()->begin()];
635 
636  solution_elem_[sbi][i_cell] = max(solution_elem_[sbi][i_cell]/n_dofs, 0.);
637  }
638  ++i_cell;
639  }
640 }
641 
642 
643 
644 
645 template<class Model>
647 {
648  //if (!Model::time_->is_current( Model::time_->marks().type_output() )) return;
649 
650 
651  START_TIMER("DG-OUTPUT");
652 
653  // gather the solution from all processors
654  data_.output_fields.set_time( this->time().step(), LimitSide::left);
655  //if (data_.output_fields.is_field_output_time(data_.output_field, this->time().step()) )
657  data_.output_fields.output(this->time().step());
658 
659  Model::output_data();
660 
661  START_TIMER("TOS-balance");
662  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
663  Model::balance_->calculate_instant(Model::subst_idx[sbi], ls[sbi]->get_solution());
664  Model::balance_->output();
665  END_TIMER("TOS-balance");
666 
667  END_TIMER("DG-OUTPUT");
668 }
669 
670 
671 template<class Model>
673 {
674  if (Model::balance_->cumulative())
675  {
676  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
677  {
678  Model::balance_->calculate_cumulative(Model::subst_idx[sbi], ls[sbi]->get_solution());
679 
680  // update source increment due to retardation
681  VecDot(ret_vec[sbi], ls[sbi]->get_solution(), &ret_sources[sbi]);
682 
683  Model::balance_->add_cumulative_source(Model::subst_idx[sbi], (ret_sources[sbi]-ret_sources_prev[sbi])/Model::time_->dt());
684  ret_sources_prev[sbi] = ret_sources[sbi];
685  }
686  }
687 }
688 
689 
690 template<class Model>
692 {
693  START_TIMER("assemble_mass");
694  Model::balance_->start_mass_assembly(Model::subst_idx);
695  assemble_mass_matrix<1>();
696  assemble_mass_matrix<2>();
697  assemble_mass_matrix<3>();
698  Model::balance_->finish_mass_assembly(Model::subst_idx);
699  END_TIMER("assemble_mass");
700 }
701 
702 
703 template<class Model> template<unsigned int dim>
705 {
706  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(), update_values | update_JxW_values | update_quadrature_points);
707  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
708  vector<LongIdx> dof_indices(ndofs);
709  PetscScalar local_mass_matrix[ndofs*ndofs], local_retardation_balance_vector[ndofs];
710  vector<PetscScalar> local_mass_balance_vector(ndofs);
711 
712  // assemble integral over elements
713  for (auto cell : feo->dh()->own_range() )
714  {
715  if (cell.dim() != dim) continue;
716  ElementAccessor<3> elm = cell.elm();
717 
718  fe_values.reinit(elm);
719  cell.get_dof_indices(dof_indices);
720 
721  Model::compute_mass_matrix_coefficient(fe_values.point_list(), elm, mm_coef);
722  Model::compute_retardation_coefficient(fe_values.point_list(), elm, ret_coef);
723 
724  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
725  {
726  // assemble the local mass matrix
727  for (unsigned int i=0; i<ndofs; i++)
728  {
729  for (unsigned int j=0; j<ndofs; j++)
730  {
731  local_mass_matrix[i*ndofs+j] = 0;
732  for (unsigned int k=0; k<qsize; k++)
733  local_mass_matrix[i*ndofs+j] += (mm_coef[k]+ret_coef[sbi][k])*fe_values.shape_value(j,k)*fe_values.shape_value(i,k)*fe_values.JxW(k);
734  }
735  }
736 
737  for (unsigned int i=0; i<ndofs; i++)
738  {
739  local_mass_balance_vector[i] = 0;
740  local_retardation_balance_vector[i] = 0;
741  for (unsigned int k=0; k<qsize; k++)
742  {
743  local_mass_balance_vector[i] += mm_coef[k]*fe_values.shape_value(i,k)*fe_values.JxW(k);
744  local_retardation_balance_vector[i] -= ret_coef[sbi][k]*fe_values.shape_value(i,k)*fe_values.JxW(k);
745  }
746  }
747 
748  Model::balance_->add_mass_matrix_values(Model::subst_idx[sbi], elm.region().bulk_idx(), dof_indices, local_mass_balance_vector);
749  ls_dt[sbi]->mat_set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), local_mass_matrix);
750  VecSetValues(ret_vec[sbi], ndofs, &(dof_indices[0]), local_retardation_balance_vector, ADD_VALUES);
751  }
752  }
753 }
754 
755 
756 
757 
758 template<class Model>
760 {
761  START_TIMER("assemble_stiffness");
762  START_TIMER("assemble_volume_integrals");
763  assemble_volume_integrals<1>();
764  assemble_volume_integrals<2>();
765  assemble_volume_integrals<3>();
766  END_TIMER("assemble_volume_integrals");
767 
768  START_TIMER("assemble_fluxes_boundary");
769  assemble_fluxes_boundary<1>();
770  assemble_fluxes_boundary<2>();
771  assemble_fluxes_boundary<3>();
772  END_TIMER("assemble_fluxes_boundary");
773 
774  START_TIMER("assemble_fluxes_elem_elem");
775  assemble_fluxes_element_element<1>();
776  assemble_fluxes_element_element<2>();
777  assemble_fluxes_element_element<3>();
778  END_TIMER("assemble_fluxes_elem_elem");
779 
780  START_TIMER("assemble_fluxes_elem_side");
781  assemble_fluxes_element_side<1>();
782  assemble_fluxes_element_side<2>();
783  assemble_fluxes_element_side<3>();
784  END_TIMER("assemble_fluxes_elem_side");
785  END_TIMER("assemble_stiffness");
786 }
787 
788 
789 
790 template<class Model>
791 template<unsigned int dim>
793 {
794  FEValues<dim,3> fv_rt(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe_rt<dim>(),
796  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
798  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
799  std::vector<LongIdx> dof_indices(ndofs);
800  vector<arma::vec3> velocity(qsize);
801  vector<vector<double> > sources_sigma(Model::n_substances(), std::vector<double>(qsize));
802  PetscScalar local_matrix[ndofs*ndofs];
803 
804  // assemble integral over elements
805  for (auto cell : feo->dh()->own_range() )
806  {
807  if (cell.dim() != dim) continue;
808  ElementAccessor<3> elm = cell.elm();
809 
810  fe_values.reinit(elm);
811  fv_rt.reinit(elm);
812  cell.get_dof_indices(dof_indices);
813 
814  calculate_velocity(elm, velocity, fv_rt);
815  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), velocity, elm, ad_coef, dif_coef);
816  Model::compute_sources_sigma(fe_values.point_list(), elm, sources_sigma);
817 
818  // assemble the local stiffness matrix
819  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
820  {
821  for (unsigned int i=0; i<ndofs; i++)
822  for (unsigned int j=0; j<ndofs; j++)
823  local_matrix[i*ndofs+j] = 0;
824 
825  for (unsigned int k=0; k<qsize; k++)
826  {
827  for (unsigned int i=0; i<ndofs; i++)
828  {
829  arma::vec3 Kt_grad_i = dif_coef[sbi][k].t()*fe_values.shape_grad(i,k);
830  double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.shape_grad(i,k));
831 
832  for (unsigned int j=0; j<ndofs; j++)
833  local_matrix[i*ndofs+j] += (arma::dot(Kt_grad_i, fe_values.shape_grad(j,k))
834  -fe_values.shape_value(j,k)*ad_dot_grad_i
835  +sources_sigma[sbi][k]*fe_values.shape_value(j,k)*fe_values.shape_value(i,k))*fe_values.JxW(k);
836  }
837  }
838  ls[sbi]->mat_set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), local_matrix);
839  }
840  }
841 }
842 
843 
844 template<class Model>
846 {
847  START_TIMER("assemble_sources");
848  Model::balance_->start_source_assembly(Model::subst_idx);
849  set_sources<1>();
850  set_sources<2>();
851  set_sources<3>();
852  Model::balance_->finish_source_assembly(Model::subst_idx);
853  END_TIMER("assemble_sources");
854 }
855 
856 template<class Model>
857 template<unsigned int dim>
859 {
860  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
862  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
863  vector<std::vector<double> > sources_conc(Model::n_substances(), std::vector<double>(qsize)),
864  sources_density(Model::n_substances(), std::vector<double>(qsize)),
865  sources_sigma(Model::n_substances(), std::vector<double>(qsize));
866  vector<LongIdx> dof_indices(ndofs);
867  PetscScalar local_rhs[ndofs];
868  vector<PetscScalar> local_source_balance_vector(ndofs), local_source_balance_rhs(ndofs);
869  double source;
870 
871  // assemble integral over elements
872  for (auto cell : feo->dh()->own_range() )
873  {
874  if (cell.dim() != dim) continue;
875  ElementAccessor<3> elm = cell.elm();
876 
877  fe_values.reinit(elm);
878  cell.get_dof_indices(dof_indices);
879 
880  Model::compute_source_coefficients(fe_values.point_list(), elm, sources_conc, sources_density, sources_sigma);
881 
882  // assemble the local stiffness matrix
883  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
884  {
885  fill_n(local_rhs, ndofs, 0);
886  local_source_balance_vector.assign(ndofs, 0);
887  local_source_balance_rhs.assign(ndofs, 0);
888 
889  // compute sources
890  for (unsigned int k=0; k<qsize; k++)
891  {
892  source = (sources_density[sbi][k] + sources_conc[sbi][k]*sources_sigma[sbi][k])*fe_values.JxW(k);
893 
894  for (unsigned int i=0; i<ndofs; i++)
895  local_rhs[i] += source*fe_values.shape_value(i,k);
896  }
897  ls[sbi]->rhs_set_values(ndofs, &(dof_indices[0]), local_rhs);
898 
899  for (unsigned int i=0; i<ndofs; i++)
900  {
901  for (unsigned int k=0; k<qsize; k++)
902  local_source_balance_vector[i] -= sources_sigma[sbi][k]*fe_values.shape_value(i,k)*fe_values.JxW(k);
903 
904  local_source_balance_rhs[i] += local_rhs[i];
905  }
906  Model::balance_->add_source_matrix_values(Model::subst_idx[sbi], elm.region().bulk_idx(), dof_indices, local_source_balance_vector);
907  Model::balance_->add_source_vec_values(Model::subst_idx[sbi], elm.region().bulk_idx(), dof_indices, local_source_balance_rhs);
908  }
909  }
910 }
911 
912 
913 
914 template<class Model>
915 template<unsigned int dim>
917 {
918  vector<FESideValues<dim,3>*> fe_values;
919  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
920  update_values);
921  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size(),
922  n_max_sides = ad_coef_edg.size();
923  vector< vector<LongIdx> > side_dof_indices;
924  PetscScalar local_matrix[ndofs*ndofs];
925  vector<vector<arma::vec3> > side_velocity(n_max_sides);
926  vector<vector<double> > dg_penalty(n_max_sides);
927  double gamma_l, omega[2], transport_flux;
928 
929  for (unsigned int sid=0; sid<n_max_sides; sid++) // Optimize: SWAP LOOPS
930  {
931  side_dof_indices.push_back( vector<LongIdx>(ndofs) );
932  fe_values.push_back(new FESideValues<dim,3>(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
934  }
935 
936  // assemble integral over sides
937  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
938  {
939  Edge *edg = &Model::mesh_->edges[feo->dh()->edge_index(iedg)];
940  if (edg->n_sides < 2 || edg->side(0)->element()->dim() != dim) continue;
941 
942  for (int sid=0; sid<edg->n_sides; sid++)
943  {
944  ElementAccessor<3> cell = Model::mesh_->element_accessor( edg->side(sid)->element().idx() );
945  feo->dh()->get_dof_indices(cell, side_dof_indices[sid]);
946  fe_values[sid]->reinit(cell, edg->side(sid)->side_idx());
947  fsv_rt.reinit(cell, edg->side(sid)->side_idx());
948  calculate_velocity(cell, side_velocity[sid], fsv_rt);
949  Model::compute_advection_diffusion_coefficients(fe_values[sid]->point_list(), side_velocity[sid], cell, ad_coef_edg[sid], dif_coef_edg[sid]);
950  dg_penalty[sid].resize(Model::n_substances());
951  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
952  dg_penalty[sid][sbi] = data_.dg_penalty[sbi].value(cell.centre(), cell);
953  }
954 
955  // fluxes and penalty
956  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
957  {
958  vector<double> fluxes(edg->n_sides);
959  for (int sid=0; sid<edg->n_sides; sid++)
960  {
961  fluxes[sid] = 0;
962  for (unsigned int k=0; k<qsize; k++)
963  fluxes[sid] += arma::dot(ad_coef_edg[sid][sbi][k], fe_values[sid]->normal_vector(k))*fe_values[sid]->JxW(k);
964  fluxes[sid] /= edg->side(sid)->measure();
965  }
966 
967  for (int s1=0; s1<edg->n_sides; s1++)
968  {
969  for (int s2=s1+1; s2<edg->n_sides; s2++)
970  {
971  OLD_ASSERT(edg->side(s1)->valid(), "Invalid side of edge.");
972  if (!feo->dh()->el_is_local( edg->side(s1)->element().idx() )
973  && !feo->dh()->el_is_local( edg->side(s2)->element().idx() )) continue;
974 
975  arma::vec3 nv = fe_values[s1]->normal_vector(0);
976 
977  // set up the parameters for DG method
978  set_DG_parameters_edge(*edg, s1, s2, qsize, dif_coef_edg[s1][sbi], dif_coef_edg[s2][sbi], fluxes, fe_values[0]->normal_vector(0), dg_penalty[s1][sbi], dg_penalty[s2][sbi], gamma_l, omega, transport_flux);
979 
980  int sd[2];
981  sd[0] = s1;
982  sd[1] = s2;
983 
984 #define AVERAGE(i,k,side_id) (fe_values[sd[side_id]]->shape_value(i,k)*0.5)
985 #define WAVERAGE(i,k,side_id) (arma::dot(dif_coef_edg[sd[side_id]][sbi][k]*fe_values[sd[side_id]]->shape_grad(i,k),nv)*omega[side_id])
986 #define JUMP(i,k,side_id) ((side_id==0?1:-1)*fe_values[sd[side_id]]->shape_value(i,k))
987 
988  // For selected pair of elements:
989  for (int n=0; n<2; n++)
990  {
991  if (!feo->dh()->el_is_local( edg->side(sd[n])->element().idx() )) continue;
992 
993  for (int m=0; m<2; m++)
994  {
995  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
996  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
997  local_matrix[i*fe_values[sd[m]]->n_dofs()+j] = 0;
998 
999  for (unsigned int k=0; k<qsize; k++)
1000  {
1001  double flux_times_JxW = transport_flux*fe_values[0]->JxW(k);
1002  double gamma_times_JxW = gamma_l*fe_values[0]->JxW(k);
1003 
1004  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
1005  {
1006  double flux_JxW_jump_i = flux_times_JxW*JUMP(i,k,n);
1007  double gamma_JxW_jump_i = gamma_times_JxW*JUMP(i,k,n);
1008  double JxW_jump_i = fe_values[0]->JxW(k)*JUMP(i,k,n);
1009  double JxW_var_wavg_i = fe_values[0]->JxW(k)*WAVERAGE(i,k,n)*dg_variant;
1010 
1011  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
1012  {
1013  int index = i*fe_values[sd[m]]->n_dofs()+j;
1014 
1015  // flux due to transport (applied on interior edges) (average times jump)
1016  local_matrix[index] += flux_JxW_jump_i*AVERAGE(j,k,m);
1017 
1018  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
1019  local_matrix[index] += gamma_JxW_jump_i*JUMP(j,k,m);
1020 
1021  // terms due to diffusion
1022  local_matrix[index] -= WAVERAGE(j,k,m)*JxW_jump_i;
1023  local_matrix[index] -= JUMP(j,k,m)*JxW_var_wavg_i;
1024  }
1025  }
1026  }
1027  ls[sbi]->mat_set_values(fe_values[sd[n]]->n_dofs(), &(side_dof_indices[sd[n]][0]), fe_values[sd[m]]->n_dofs(), &(side_dof_indices[sd[m]][0]), local_matrix);
1028  }
1029  }
1030 #undef AVERAGE
1031 #undef WAVERAGE
1032 #undef JUMP
1033  }
1034  }
1035  }
1036  }
1037 
1038  for (unsigned int i=0; i<n_max_sides; i++)
1039  {
1040  delete fe_values[i];
1041  }
1042 }
1043 
1044 
1045 template<class Model>
1046 template<unsigned int dim>
1048 {
1049  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1051  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1052  update_values);
1053  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1054  std::vector<LongIdx> side_dof_indices(ndofs);
1055  PetscScalar local_matrix[ndofs*ndofs];
1056  vector<arma::vec3> side_velocity;
1057  vector<double> robin_sigma(qsize);
1058  vector<double> csection(qsize);
1059  double gamma_l;
1060 
1061  // assemble boundary integral
1062  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1063  {
1064  Edge *edg = &Model::mesh_->edges[feo->dh()->edge_index(iedg)];
1065  if (edg->n_sides > 1) continue;
1066  // check spatial dimension
1067  if (edg->side(0)->dim() != dim-1) continue;
1068  // skip edges lying not on the boundary
1069  if (edg->side(0)->cond() == NULL) continue;
1070 
1071  SideIter side = edg->side(0);
1072  auto cell = feo->dh()->cell_accessor_from_element( side->element().idx() );
1073  ElementAccessor<3> elm_acc = cell.elm();
1074  cell.get_dof_indices(side_dof_indices);
1075  fe_values_side.reinit(elm_acc, side->side_idx());
1076  fsv_rt.reinit(elm_acc, side->side_idx());
1077 
1078  calculate_velocity(elm_acc, side_velocity, fsv_rt);
1079  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), side_velocity, elm_acc, ad_coef, dif_coef);
1080  arma::uvec bc_type;
1081  Model::get_bc_type(side->cond()->element_accessor(), bc_type);
1082  data_.cross_section.value_list(fe_values_side.point_list(), elm_acc, csection);
1083 
1084  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1085  {
1086  for (unsigned int i=0; i<ndofs; i++)
1087  for (unsigned int j=0; j<ndofs; j++)
1088  local_matrix[i*ndofs+j] = 0;
1089 
1090  // On Neumann boundaries we have only term from integrating by parts the advective term,
1091  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
1092  double side_flux = 0;
1093  for (unsigned int k=0; k<qsize; k++)
1094  side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
1095  double transport_flux = side_flux/side->measure();
1096 
1097  if (bc_type[sbi] == AdvectionDiffusionModel::abc_dirichlet)
1098  {
1099  // set up the parameters for DG method
1100  set_DG_parameters_boundary(side, qsize, dif_coef[sbi], transport_flux, fe_values_side.normal_vector(0), data_.dg_penalty[sbi].value(elm_acc.centre(), elm_acc), gamma_l);
1101  gamma[sbi][side->cond_idx()] = gamma_l;
1102  transport_flux += gamma_l;
1103  }
1104 
1105  // fluxes and penalty
1106  for (unsigned int k=0; k<qsize; k++)
1107  {
1108  double flux_times_JxW;
1109  if (bc_type[sbi] == AdvectionDiffusionModel::abc_total_flux)
1110  {
1111  Model::get_flux_bc_sigma(sbi, fe_values_side.point_list(), side->cond()->element_accessor(), robin_sigma);
1112  flux_times_JxW = csection[k]*robin_sigma[k]*fe_values_side.JxW(k);
1113  }
1114  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_diffusive_flux)
1115  {
1116  Model::get_flux_bc_sigma(sbi, fe_values_side.point_list(), side->cond()->element_accessor(), robin_sigma);
1117  flux_times_JxW = (transport_flux + csection[k]*robin_sigma[k])*fe_values_side.JxW(k);
1118  }
1119  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
1120  flux_times_JxW = 0;
1121  else
1122  flux_times_JxW = transport_flux*fe_values_side.JxW(k);
1123 
1124  for (unsigned int i=0; i<ndofs; i++)
1125  {
1126  for (unsigned int j=0; j<ndofs; j++)
1127  {
1128  // flux due to advection and penalty
1129  local_matrix[i*ndofs+j] += flux_times_JxW*fe_values_side.shape_value(i,k)*fe_values_side.shape_value(j,k);
1130 
1131  // flux due to diffusion (only on dirichlet and inflow boundary)
1132  if (bc_type[sbi] == AdvectionDiffusionModel::abc_dirichlet)
1133  local_matrix[i*ndofs+j] -= (arma::dot(dif_coef[sbi][k]*fe_values_side.shape_grad(j,k),fe_values_side.normal_vector(k))*fe_values_side.shape_value(i,k)
1134  + arma::dot(dif_coef[sbi][k]*fe_values_side.shape_grad(i,k),fe_values_side.normal_vector(k))*fe_values_side.shape_value(j,k)*dg_variant
1135  )*fe_values_side.JxW(k);
1136  }
1137  }
1138  }
1139 
1140  ls[sbi]->mat_set_values(ndofs, &(side_dof_indices[0]), ndofs, &(side_dof_indices[0]), local_matrix);
1141  }
1142  }
1143 }
1144 
1145 
1146 template<class Model>
1147 template<unsigned int dim>
1149 {
1150 
1151  if (dim == 1) return;
1152  FEValues<dim-1,3> fe_values_vb(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe<dim-1>(),
1154  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1156  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1157  update_values);
1158  FEValues<dim-1,3> fv_rt(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe_rt<dim-1>(),
1159  update_values);
1160 
1161  vector<FEValuesSpaceBase<3>*> fv_sb(2);
1162  const unsigned int ndofs = feo->fe<dim>()->n_dofs(); // number of local dofs
1163  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
1164  int side_dof_indices[2*ndofs];
1165  vector<LongIdx> indices(ndofs);
1166  unsigned int n_dofs[2], n_indices;
1167  vector<arma::vec3> velocity_higher, velocity_lower;
1168  vector<double> frac_sigma(qsize);
1169  vector<double> csection_lower(qsize), csection_higher(qsize);
1170  PetscScalar local_matrix[4*ndofs*ndofs];
1171  double comm_flux[2][2];
1172 
1173  // index 0 = element with lower dimension,
1174  // index 1 = side of element with higher dimension
1175  fv_sb[0] = &fe_values_vb;
1176  fv_sb[1] = &fe_values_side;
1177 
1178  // assemble integral over sides
1179  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
1180  {
1181  Neighbour *nb = &Model::mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
1182  // skip neighbours of different dimension
1183  if (nb->element()->dim() != dim-1) continue;
1184 
1185  ElementAccessor<3> cell_sub = Model::mesh_->element_accessor( nb->element().idx() );
1186  n_indices = feo->dh()->get_dof_indices(cell_sub, indices);
1187  for(unsigned int i=0; i<n_indices; ++i) {
1188  side_dof_indices[i] = indices[i];
1189  }
1190  fe_values_vb.reinit(cell_sub);
1191  n_dofs[0] = fv_sb[0]->n_dofs();
1192 
1193  ElementAccessor<3> cell = Model::mesh_->element_accessor( nb->side()->element().idx() );
1194  n_indices = feo->dh()->get_dof_indices(cell, indices);
1195  for(unsigned int i=0; i<n_indices; ++i) {
1196  side_dof_indices[i+n_dofs[0]] = indices[i];
1197  }
1198  fe_values_side.reinit(cell, nb->side()->side_idx());
1199  n_dofs[1] = fv_sb[1]->n_dofs();
1200 
1201  // Element id's for testing if they belong to local partition.
1202  int element_id[2];
1203  element_id[0] = cell_sub.idx();
1204  element_id[1] = cell.idx();
1205 
1206  fsv_rt.reinit(cell, nb->side()->side_idx());
1207  fv_rt.reinit(cell_sub);
1208  calculate_velocity(cell, velocity_higher, fsv_rt);
1209  calculate_velocity(cell_sub, velocity_lower, fv_rt);
1210  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_lower, cell_sub, ad_coef_edg[0], dif_coef_edg[0]);
1211  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_higher, cell, ad_coef_edg[1], dif_coef_edg[1]);
1212  data_.cross_section.value_list(fe_values_vb.point_list(), cell_sub, csection_lower);
1213  data_.cross_section.value_list(fe_values_vb.point_list(), cell, csection_higher);
1214 
1215  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++) // Optimize: SWAP LOOPS
1216  {
1217  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
1218  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
1219  local_matrix[i*(n_dofs[0]+n_dofs[1])+j] = 0;
1220 
1221  data_.fracture_sigma[sbi].value_list(fe_values_vb.point_list(), cell_sub, frac_sigma);
1222 
1223  // set transmission conditions
1224  for (unsigned int k=0; k<qsize; k++)
1225  {
1226  /* The communication flux has two parts:
1227  * - "diffusive" term containing sigma
1228  * - "advective" term representing usual upwind
1229  *
1230  * The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
1231  * than b and A in the manual.
1232  * In calculation of sigma there appears one more csection_lower in the denominator.
1233  */
1234  double sigma = frac_sigma[k]*arma::dot(dif_coef_edg[0][sbi][k]*fe_values_side.normal_vector(k),fe_values_side.normal_vector(k))*
1235  2*csection_higher[k]*csection_higher[k]/(csection_lower[k]*csection_lower[k]);
1236 
1237  double transport_flux = arma::dot(ad_coef_edg[1][sbi][k], fe_values_side.normal_vector(k));
1238 
1239  comm_flux[0][0] = (sigma-min(0.,transport_flux))*fv_sb[0]->JxW(k);
1240  comm_flux[0][1] = -(sigma-min(0.,transport_flux))*fv_sb[0]->JxW(k);
1241  comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1242  comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1243 
1244  for (int n=0; n<2; n++)
1245  {
1246  if (!feo->dh()->el_is_local(element_id[n])) continue;
1247 
1248  for (unsigned int i=0; i<n_dofs[n]; i++)
1249  for (int m=0; m<2; m++)
1250  for (unsigned int j=0; j<n_dofs[m]; j++)
1251  local_matrix[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
1252  comm_flux[m][n]*fv_sb[m]->shape_value(j,k)*fv_sb[n]->shape_value(i,k);
1253  }
1254  }
1255  ls[sbi]->mat_set_values(n_dofs[0]+n_dofs[1], side_dof_indices, n_dofs[0]+n_dofs[1], side_dof_indices, local_matrix);
1256  }
1257  }
1258 
1259 }
1260 
1261 
1262 
1263 
1264 
1265 
1266 template<class Model>
1268 {
1269  START_TIMER("assemble_bc");
1270  Model::balance_->start_flux_assembly(Model::subst_idx);
1271  set_boundary_conditions<1>();
1272  set_boundary_conditions<2>();
1273  set_boundary_conditions<3>();
1274  Model::balance_->finish_flux_assembly(Model::subst_idx);
1275  END_TIMER("assemble_bc");
1276 }
1277 
1278 
1279 template<class Model>
1280 template<unsigned int dim>
1282 {
1283  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1285  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1286  update_values);
1287  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1288  vector<LongIdx> side_dof_indices(ndofs);
1289  unsigned int loc_b=0;
1290  double local_rhs[ndofs];
1291  vector<PetscScalar> local_flux_balance_vector(ndofs);
1292  PetscScalar local_flux_balance_rhs;
1293  vector<double> bc_values(qsize);
1294  vector<double> bc_fluxes(qsize),
1295  bc_sigma(qsize),
1296  bc_ref_values(qsize);
1297  vector<double> csection(qsize);
1298  vector<arma::vec3> velocity;
1299 
1300  for (auto cell : feo->dh()->own_range() )
1301  {
1302  if (cell.elm()->boundary_idx_ == nullptr) continue;
1303 
1304  for (unsigned int si=0; si<cell.elm()->n_sides(); si++)
1305  {
1306  const Edge *edg = cell.elm().side(si)->edge();
1307  if (edg->n_sides > 1) continue;
1308  // skip edges lying not on the boundary
1309  if (edg->side(0)->cond() == NULL) continue;
1310 
1311  if (edg->side(0)->dim() != dim-1)
1312  {
1313  if (edg->side(0)->cond() != nullptr) ++loc_b;
1314  continue;
1315  }
1316 
1317  SideIter side = edg->side(0);
1318  ElementAccessor<3> elm = Model::mesh_->element_accessor( side->element().idx() );
1319  ElementAccessor<3> ele_acc = side->cond()->element_accessor();
1320 
1321  arma::uvec bc_type;
1322  Model::get_bc_type(ele_acc, bc_type);
1323 
1324  fe_values_side.reinit(elm, side->side_idx());
1325  fsv_rt.reinit(elm, side->side_idx());
1326  calculate_velocity(elm, velocity, fsv_rt);
1327 
1328  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), velocity, side->element(), ad_coef, dif_coef);
1329  data_.cross_section.value_list(fe_values_side.point_list(), side->element(), csection);
1330 
1331  auto dh_cell = feo->dh()->cell_accessor_from_element( side->element().idx() );
1332  dh_cell.get_dof_indices(side_dof_indices);
1333 
1334  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1335  {
1336  fill_n(local_rhs, ndofs, 0);
1337  local_flux_balance_vector.assign(ndofs, 0);
1338  local_flux_balance_rhs = 0;
1339 
1340  // The b.c. data are fetched for all possible b.c. types since we allow
1341  // different bc_type for each substance.
1342  data_.bc_dirichlet_value[sbi].value_list(fe_values_side.point_list(), ele_acc, bc_values);
1343 
1344  double side_flux = 0;
1345  for (unsigned int k=0; k<qsize; k++)
1346  side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
1347  double transport_flux = side_flux/side->measure();
1348 
1349  if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
1350  {
1351  for (unsigned int k=0; k<qsize; k++)
1352  {
1353  double bc_term = -transport_flux*bc_values[k]*fe_values_side.JxW(k);
1354  for (unsigned int i=0; i<ndofs; i++)
1355  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k);
1356  }
1357  for (unsigned int i=0; i<ndofs; i++)
1358  local_flux_balance_rhs -= local_rhs[i];
1359  }
1360  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_dirichlet)
1361  {
1362  for (unsigned int k=0; k<qsize; k++)
1363  {
1364  double bc_term = gamma[sbi][side->cond_idx()]*bc_values[k]*fe_values_side.JxW(k);
1365  arma::vec3 bc_grad = -bc_values[k]*fe_values_side.JxW(k)*dg_variant*(arma::trans(dif_coef[sbi][k])*fe_values_side.normal_vector(k));
1366  for (unsigned int i=0; i<ndofs; i++)
1367  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k)
1368  + arma::dot(bc_grad,fe_values_side.shape_grad(i,k));
1369  }
1370  for (unsigned int k=0; k<qsize; k++)
1371  {
1372  for (unsigned int i=0; i<ndofs; i++)
1373  {
1374  local_flux_balance_vector[i] += (arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.shape_value(i,k)
1375  - arma::dot(dif_coef[sbi][k]*fe_values_side.shape_grad(i,k),fe_values_side.normal_vector(k))
1376  + gamma[sbi][side->cond_idx()]*fe_values_side.shape_value(i,k))*fe_values_side.JxW(k);
1377  }
1378  }
1379  if (Model::time_->tlevel() > 0)
1380  for (unsigned int i=0; i<ndofs; i++)
1381  local_flux_balance_rhs -= local_rhs[i];
1382  }
1383  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_total_flux)
1384  {
1385  Model::get_flux_bc_data(sbi, fe_values_side.point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1386  for (unsigned int k=0; k<qsize; k++)
1387  {
1388  double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.JxW(k);
1389  for (unsigned int i=0; i<ndofs; i++)
1390  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k);
1391  }
1392 
1393  for (unsigned int i=0; i<ndofs; i++)
1394  {
1395  for (unsigned int k=0; k<qsize; k++)
1396  local_flux_balance_vector[i] += csection[k]*bc_sigma[k]*fe_values_side.JxW(k)*fe_values_side.shape_value(i,k);
1397  local_flux_balance_rhs -= local_rhs[i];
1398  }
1399  }
1400  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_diffusive_flux)
1401  {
1402  Model::get_flux_bc_data(sbi, fe_values_side.point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1403  for (unsigned int k=0; k<qsize; k++)
1404  {
1405  double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.JxW(k);
1406  for (unsigned int i=0; i<ndofs; i++)
1407  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k);
1408  }
1409 
1410  for (unsigned int i=0; i<ndofs; i++)
1411  {
1412  for (unsigned int k=0; k<qsize; k++)
1413  local_flux_balance_vector[i] += csection[k]*(arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k)) + bc_sigma[k])*fe_values_side.JxW(k)*fe_values_side.shape_value(i,k);
1414  local_flux_balance_rhs -= local_rhs[i];
1415  }
1416  }
1417  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux >= 0)
1418  {
1419  for (unsigned int k=0; k<qsize; k++)
1420  {
1421  for (unsigned int i=0; i<ndofs; i++)
1422  local_flux_balance_vector[i] += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k)*fe_values_side.shape_value(i,k);
1423  }
1424  }
1425  ls[sbi]->rhs_set_values(ndofs, &(side_dof_indices[0]), local_rhs);
1426 
1427  Model::balance_->add_flux_matrix_values(Model::subst_idx[sbi], loc_b, side_dof_indices, local_flux_balance_vector);
1428  Model::balance_->add_flux_vec_value(Model::subst_idx[sbi], loc_b, local_flux_balance_rhs);
1429  }
1430  ++loc_b;
1431  }
1432  }
1433 }
1434 
1435 
1436 
1437 template<class Model>
1438 template<unsigned int dim>
1440  vector<arma::vec3> &velocity,
1441  FEValuesBase<dim,3> &fv)
1442 {
1443  OLD_ASSERT(cell->dim() == dim, "Element dimension mismatch!");
1444 
1445  velocity.resize(fv.n_points());
1446 
1447  for (unsigned int k=0; k<fv.n_points(); k++)
1448  {
1449  velocity[k].zeros();
1450  for (unsigned int sid=0; sid<cell->n_sides(); sid++)
1451  for (unsigned int c=0; c<3; ++c)
1452  velocity[k][c] += fv.shape_value_component(sid,k,c) * Model::mh_dh->side_flux( *(cell.side(sid)) );
1453  }
1454 }
1455 
1456 
1457 // return the ratio of longest and shortest edge
1459 {
1460  double h_max = 0, h_min = numeric_limits<double>::infinity();
1461  for (unsigned int i=0; i<e->n_nodes(); i++)
1462  for (unsigned int j=i+1; j<e->n_nodes(); j++)
1463  {
1464  h_max = max(h_max, e.node(i)->distance(*e.node(j)));
1465  h_min = min(h_min, e.node(i)->distance(*e.node(j)));
1466  }
1467  return h_max/h_min;
1468 }
1469 
1470 
1471 
1472 template<class Model>
1474  const int s1,
1475  const int s2,
1476  const int K_size,
1477  const vector<arma::mat33> &K1,
1478  const vector<arma::mat33> &K2,
1479  const vector<double> &fluxes,
1480  const arma::vec3 &normal_vector,
1481  const double alpha1,
1482  const double alpha2,
1483  double &gamma,
1484  double *omega,
1485  double &transport_flux)
1486 {
1487  double delta[2];
1488  double h = 0;
1489  double local_alpha = 0;
1490 
1491  OLD_ASSERT(edg.side(s1)->valid(), "Invalid side of an edge.");
1492  SideIter s = edg.side(s1);
1493 
1494  // calculate the side diameter
1495  if (s->dim() == 0)
1496  {
1497  h = 1;
1498  }
1499  else
1500  {
1501  for (unsigned int i=0; i<s->n_nodes(); i++)
1502  for (unsigned int j=i+1; j<s->n_nodes(); j++)
1503  h = max(h, s->node(i)->distance(*s->node(j).node()));
1504  }
1505 
1506  double aniso1 = elem_anisotropy(edg.side(s1)->element());
1507  double aniso2 = elem_anisotropy(edg.side(s2)->element());
1508 
1509 
1510  // calculate the total in- and out-flux through the edge
1511  double pflux = 0, nflux = 0;
1512  for (int i=0; i<edg.n_sides; i++)
1513  {
1514  if (fluxes[i] > 0)
1515  pflux += fluxes[i];
1516  else
1517  nflux += fluxes[i];
1518  }
1519 
1520  // calculate the flux from s1 to s2
1521  if (fluxes[s2] > 0 && fluxes[s1] < 0 && s1 < s2)
1522  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
1523  else if (fluxes[s2] < 0 && fluxes[s1] > 0 && s1 < s2)
1524  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
1525  else if (s1==s2)
1526  transport_flux = fluxes[s1];
1527  else
1528  transport_flux = 0;
1529 
1530  gamma = 0.5*fabs(transport_flux);
1531 
1532 
1533  // determine local DG penalty
1534  local_alpha = max(alpha1, alpha2);
1535 
1536  if (s1 == s2)
1537  {
1538  omega[0] = 1;
1539 
1540  // delta is set to the average value of Kn.n on the side
1541  delta[0] = 0;
1542  for (int k=0; k<K_size; k++)
1543  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1544  delta[0] /= K_size;
1545 
1546  gamma += local_alpha/h*aniso1*aniso2*delta[0];
1547  }
1548  else
1549  {
1550  delta[0] = 0;
1551  delta[1] = 0;
1552  for (int k=0; k<K_size; k++)
1553  {
1554  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1555  delta[1] += dot(K2[k]*normal_vector,normal_vector);
1556  }
1557  delta[0] /= K_size;
1558  delta[1] /= K_size;
1559 
1560  double delta_sum = delta[0] + delta[1];
1561 
1562 // if (delta_sum > numeric_limits<double>::epsilon())
1563  if (fabs(delta_sum) > 0)
1564  {
1565  omega[0] = delta[1]/delta_sum;
1566  omega[1] = delta[0]/delta_sum;
1567  gamma += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
1568  }
1569  else
1570  for (int i=0; i<2; i++) omega[i] = 0;
1571  }
1572 }
1573 
1574 
1575 
1576 
1577 
1578 
1579 template<class Model>
1581  const int K_size,
1582  const vector<arma::mat33> &K,
1583  const double flux,
1584  const arma::vec3 &normal_vector,
1585  const double alpha,
1586  double &gamma)
1587 {
1588  double delta = 0, h = 0;
1589 
1590  // calculate the side diameter
1591  if (side->dim() == 0)
1592  {
1593  h = 1;
1594  }
1595  else
1596  {
1597  for (unsigned int i=0; i<side->n_nodes(); i++)
1598  for (unsigned int j=i+1; j<side->n_nodes(); j++)
1599  h = max(h, side->node(i)->distance( *side->node(j).node() ));
1600  }
1601 
1602  // delta is set to the average value of Kn.n on the side
1603  for (int k=0; k<K_size; k++)
1604  delta += dot(K[k]*normal_vector,normal_vector);
1605  delta /= K_size;
1606 
1607  gamma = 0.5*fabs(flux) + alpha/h*delta*elem_anisotropy(side->element());
1608 }
1609 
1610 
1611 
1612 
1613 
1614 template<class Model>
1616 {
1617  START_TIMER("set_init_cond");
1618  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1619  ls[sbi]->start_allocation();
1620  prepare_initial_condition<1>();
1621  prepare_initial_condition<2>();
1622  prepare_initial_condition<3>();
1623 
1624  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1625  ls[sbi]->start_add_assembly();
1626  prepare_initial_condition<1>();
1627  prepare_initial_condition<2>();
1628  prepare_initial_condition<3>();
1629 
1630  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1631  {
1632  ls[sbi]->finish_assembly();
1633  ls[sbi]->solve();
1634  }
1635  END_TIMER("set_init_cond");
1636 }
1637 
1638 template<class Model>
1639 template<unsigned int dim>
1641 {
1642  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1644  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1645  std::vector<LongIdx> dof_indices(ndofs);
1646  double matrix[ndofs*ndofs], rhs[ndofs];
1647  std::vector<std::vector<double> > init_values(Model::n_substances());
1648 
1649  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++) // Optimize: SWAP LOOPS
1650  init_values[sbi].resize(qsize);
1651 
1652  for (auto cell : feo->dh()->own_range() )
1653  {
1654  if (cell.dim() != dim) continue;
1655  ElementAccessor<3> elem = cell.elm();
1656 
1657  cell.get_dof_indices(dof_indices);
1658  fe_values.reinit(elem);
1659 
1660  Model::compute_init_cond(fe_values.point_list(), elem, init_values);
1661 
1662  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1663  {
1664  for (unsigned int i=0; i<ndofs; i++)
1665  {
1666  rhs[i] = 0;
1667  for (unsigned int j=0; j<ndofs; j++)
1668  matrix[i*ndofs+j] = 0;
1669  }
1670 
1671  for (unsigned int k=0; k<qsize; k++)
1672  {
1673  double rhs_term = init_values[sbi][k]*fe_values.JxW(k);
1674 
1675  for (unsigned int i=0; i<ndofs; i++)
1676  {
1677  for (unsigned int j=0; j<ndofs; j++)
1678  matrix[i*ndofs+j] += fe_values.shape_value(i,k)*fe_values.shape_value(j,k)*fe_values.JxW(k);
1679 
1680  rhs[i] += fe_values.shape_value(i,k)*rhs_term;
1681  }
1682  }
1683  ls[sbi]->set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), matrix, rhs);
1684  }
1685  }
1686 }
1687 
1688 
1689 template<class Model>
1691 {
1692  el_4_loc = Model::mesh_->get_el_4_loc();
1693  el_ds = Model::mesh_->get_el_ds();
1694 }
1695 
1696 
1697 template<class Model>
1699 {
1700  if (solution_changed)
1701  {
1702  unsigned int i_cell=0;
1703  for (auto cell : feo->dh()->own_range() )
1704  {
1705 
1706  unsigned int n_dofs;
1707  switch (cell.dim())
1708  {
1709  case 1:
1710  n_dofs = feo->fe<1>()->n_dofs();
1711  break;
1712  case 2:
1713  n_dofs = feo->fe<2>()->n_dofs();
1714  break;
1715  case 3:
1716  n_dofs = feo->fe<3>()->n_dofs();
1717  break;
1718  }
1719 
1720  std::vector<LongIdx> dof_indices(n_dofs);
1721  cell.get_dof_indices(dof_indices);
1722 
1723  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1724  {
1725  double old_average = 0;
1726  for (unsigned int j=0; j<n_dofs; ++j)
1727  old_average += ls[sbi]->get_solution_array()[dof_indices[j]-feo->dh()->distr()->begin()];
1728  old_average /= n_dofs;
1729 
1730  for (unsigned int j=0; j<n_dofs; ++j)
1731  ls[sbi]->get_solution_array()[dof_indices[j]-feo->dh()->distr()->begin()] += solution_elem_[sbi][i_cell] - old_average;
1732  }
1733  ++i_cell;
1734  }
1735  }
1736  // update mass_vec for the case that mass matrix changes in next time step
1737  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1738  MatMult(*(ls_dt[sbi]->get_matrix()), ls[sbi]->get_solution(), mass_vec[sbi]);
1739 }
1740 
1741 template<class Model>
1743 {
1744  return Model::mesh_->get_row_4_el();
1745 }
1746 
1747 
1748 
1749 
1750 
1751 
1752 
1753 
1755 template class TransportDG<HeatTransferModel>;
1756 
1757 
1758 
1759 
TimeGovernor & time()
Definition: equation.hh:148
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
const Edge * edge() const
Definition: side_impl.hh:66
Input::Record input_rec
Record with input specification.
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:87
FieldSet * eq_data_
Definition: equation.hh:232
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.hh:336
static auto subdomain(Mesh &mesh) -> IndexField
void assemble_fluxes_element_element()
Assembles the fluxes between elements of the same dimension.
void set_sources()
Assembles the right hand side due to volume sources.
FiniteElement< dim > * fe()
void set_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
void get_par_info(LongIdx *&el_4_loc, Distribution *&el_ds)
unsigned int n_nodes() const
Definition: elements.h:129
Transformed quadrature weight for cell sides.
vector< double > mm_coef
Mass matrix coefficients.
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
void assemble_fluxes_element_side()
Assembles the fluxes between elements of different dimensions.
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:49
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:604
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
unsigned int side_idx() const
Definition: side_impl.hh:82
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)
double distance(const Node &n2) const
Definition: nodes.hh:85
void set_from_input(const Input::Record in_rec) override
virtual void start_add_assembly()
Definition: linsys.hh:322
void assemble_mass_matrix()
Assembles the mass matrix.
void output(TimeStep step)
Boundary * cond() const
Definition: side_impl.hh:71
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:264
vector< VectorMPI > output_vec
Array for storing the output solution data.
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
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
int dg_variant
DG variant ((non-)symmetric/incomplete.
void prepare_initial_condition()
Assembles the auxiliary linear system to calculate the initial solution as L^2-projection of the pres...
Definition: mesh.h:80
Fields computed from the mesh data.
virtual void start_allocation()
Definition: linsys.hh:314
void set_initial_condition()
Sets the initial condition.
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:147
virtual void finish_assembly()=0
int n_sides
Definition: edges.h:36
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
Definition: edges.h:26
bool valid() const
Definition: side_impl.hh:92
SideIter side(const unsigned int loc_index)
Definition: accessors.hh:137
#define AVERAGE(i, k, side_id)
std::shared_ptr< DOFHandlerMultiDim > dh()
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
vector< vector< vector< arma::mat33 > > > dif_coef_edg
Diffusion coefficients on edges.
Discontinuous Galerkin method for equation of transport with dispersion.
Class for declaration of the integral input data.
Definition: type_base.hh:490
const Vec & get_solution(unsigned int sbi)
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
ElementAccessor< 3 > element() const
Definition: side_impl.hh:53
vector< double > ret_sources_prev
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
static const Input::Type::Selection & get_dg_variant_selection_input_type()
Input type for the DG variant selection.
Definition: transport_dg.cc:53
unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.hh:407
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()
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
static constexpr bool value
Definition: json.hpp:87
Symmetric Gauss-Legendre quadrature formulae on simplices.
#define WAVERAGE(i, k, side_id)
TransportDG(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
vector< vector< double > > ret_coef
Retardation coefficient due to sorption.
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors.hh:285
unsigned int dim() const
Definition: elements.h:124
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:368
#define OLD_ASSERT(...)
Definition: global_defs.h:131
ElementAccessor< 3 > element()
Definition: neighbours.h:163
LongIdx * get_row_4_el()
void initialize() override
vector< vector< vector< arma::vec3 > > > ad_coef_edg
Advection coefficients on edges.
Transformed quadrature points.
void output_data()
Postprocesses the solution and writes to output file.
virtual PetscErrorCode set_rhs(Vec &rhs)
Definition: linsys.hh:255
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
void preallocate()
unsigned int dim() const
Definition: side_impl.hh:38
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
FEObjects * feo
Finite element objects.
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)
const Vec & get_solution()
Definition: linsys.hh:282
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:60
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
SideIter side()
Definition: neighbours.h:147
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
unsigned int n_sides() const
Definition: elements.h:135
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.
Mesh * mesh_
Definition: equation.hh:223
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:501
Field< 3, FieldValue< 3 >::Scalar > region_id
double elem_anisotropy(ElementAccessor< 3 > e)
unsigned int dg_order
Polynomial order of finite elements.
void output_vector_gather()
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
Shape function gradients.
Definition: update_flags.hh:95
FLOW123D_FORCE_LINK_IN_CHILD(concentrationTransportModel)
void set_solution(double *sol_array)
Definition: linsys.hh:290
NodeAccessor< 3 > node(unsigned int i) const
Definition: side_impl.hh:47
double measure() const
Definition: sides.cc:30
double shape_value(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:224
Region region() const
Definition: accessors.hh:95
Normal vectors.
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:273
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).
#define MPI_Comm_rank
Definition: mpi.h:236
std::vector< Mat > stiffness_matrix
The stiffness matrix.
FieldCommon & description(const string &description)
void set_values(int nrow, int *rows, int ncol, int *cols, PetscScalar *mat_vals, PetscScalar *rhs_vals)
Set values in the system matrix and values in the right-hand side vector on corresponding rows...
Definition: linsys.hh:367
EqData data_
Field data for model parameters.
bool allocation_done
Indicates whether matrices have been preallocated.
unsigned int cond_idx() const
Definition: side_impl.hh:76
std::vector< std::vector< double > > gamma
Penalty parameters.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: transport_dg.cc:79
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
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)
MappingP1< dim, 3 > * mapping()
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:242
vector< double > ret_sources
Temporary values of increments due to retardation (e.g. sorption)
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no)
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:233
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
FieldCommon & name(const string &name)
virtual SolveInfo solve()=0
#define END_TIMER(tag)
Ends a timer with specified tag.
Discontinuous Galerkin method for equation of transport with dispersion.
void calculate_velocity(const ElementAccessor< 3 > &cell, std::vector< arma::vec3 > &velocity, FEValuesBase< dim, 3 > &fv)
Calculates the velocity field on a given dim dimensional cell.
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.hh:357
Quadrature< dim > * q()
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:32
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
Definition: linsys.hh:246
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
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Base class for FEValues and FESideValues.
Definition: fe_values.hh:37
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
~TransportDG()
Destructor.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:111
LinSys ** ls
Linear algebra system for the transport equation.
FiniteElement< dim > * fe_rt()
void set_DG_parameters_edge(const Edge &edg, const int s1, const int s2, const int K_size, const std::vector< arma::mat33 > &K1, const std::vector< arma::mat33 > &K2, const std::vector< double > &fluxes, const arma::vec3 &normal_vector, const double alpha1, const double alpha2, double &gamma, double *omega, double &transport_flux)
Calculates the dispersivity (diffusivity) tensor from the velocity field.
SideIter side(const unsigned int i) const
Definition: edges.h:31
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
#define JUMP(i, k, side_id)
Template for classes storing finite set of named values.
void set_DG_parameters_boundary(const SideIter 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.
Definitions of Raviart-Thomas finite elements.
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:48
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)
Definition: field.hh:56
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:528
unsigned int n_nodes() const
Definition: side_impl.hh:34
Transformed quadrature weights.
Calculates finite element data on a side.
Definition: fe_values.hh:577
const Node * node(unsigned int ni) const
Definition: accessors.hh:145