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