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