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