Flow123d  last_with_con_2.0.0-663-gd0e2296
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"
38 #include "io/equation_output.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 
230  // add all input fields to the output list
231  output_fields += *this;
232 
233 }
234 
235 template<class Model>
237  : Model(init_mesh, in_rec),
238  input_rec(in_rec),
239  allocation_done(false)
240 {
241  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
242  // due to constexpr optimization.
243  START_TIMER(Model::ModelEqData::name());
244  // Check that Model is derived from AdvectionDiffusionModel.
246 
247  this->eq_data_ = &data_;
248 
249 
250  // Set up physical parameters.
251  data_.set_mesh(init_mesh);
252  data_.region_id = GenericField<3>::region_id(*Model::mesh_);
253 
254 
255  // DG variant and order
256  dg_variant = in_rec.val<DGVariant>("dg_variant");
257  dg_order = in_rec.val<unsigned int>("dg_order");
258 
259  Model::init_from_input(in_rec);
260 
261  // create finite element structures and distribute DOFs
262  feo = new FEObjects(Model::mesh_, dg_order);
263  //DebugOut().fmt("TDG: solution size {}\n", feo->dh()->n_global_dofs());
264 
265 }
266 
267 
268 template<class Model>
270 {
271  data_.set_components(Model::substances_.names());
272  data_.set_input_list( input_rec.val<Input::Array>("input_fields") );
273 
274  // DG stabilization parameters on boundary edges
275  gamma.resize(Model::n_substances());
276  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
277  gamma[sbi].resize(Model::mesh_->boundary_.size());
278 
279  // Resize coefficient arrays
280  int qsize = max(feo->q<0>()->size(), max(feo->q<1>()->size(), max(feo->q<2>()->size(), feo->q<3>()->size())));
281  int max_edg_sides = max(Model::mesh_->max_edge_sides(1), max(Model::mesh_->max_edge_sides(2), Model::mesh_->max_edge_sides(3)));
282  mm_coef.resize(qsize);
283  ret_coef.resize(Model::n_substances());
284  sorption_sources.resize(Model::n_substances());
285  ad_coef.resize(Model::n_substances());
286  dif_coef.resize(Model::n_substances());
287  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
288  {
289  ret_coef[sbi].resize(qsize);
290  ad_coef[sbi].resize(qsize);
291  dif_coef[sbi].resize(qsize);
292  }
293  ad_coef_edg.resize(max_edg_sides);
294  dif_coef_edg.resize(max_edg_sides);
295  for (int sd=0; sd<max_edg_sides; sd++)
296  {
297  ad_coef_edg[sd].resize(Model::n_substances());
298  dif_coef_edg[sd].resize(Model::n_substances());
299  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
300  {
301  ad_coef_edg[sd][sbi].resize(qsize);
302  dif_coef_edg[sd][sbi].resize(qsize);
303  }
304  }
305 
306  output_vec.resize(Model::n_substances());
307  //output_solution.resize(Model::n_substances());
308  int rank;
309  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
310  unsigned int output_vector_size= (rank==0)?feo->dh()->n_global_dofs():0;
311  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
312  {
313  // for each substance we allocate output array and vector
314  //output_solution[sbi] = new double[feo->dh()->n_global_dofs()];
315  VecCreateSeq(PETSC_COMM_SELF, output_vector_size, &output_vec[sbi]);
316  }
317  data_.output_field.set_components(Model::substances_.names());
318  data_.output_field.set_mesh(*Model::mesh_);
319  data_.output_type(OutputTime::CORNER_DATA);
320 
321  data_.output_field.setup_components();
322  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
323  {
324  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
325  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar> > output_field_ptr(new FieldFE<3, FieldValue<3>::Scalar>);
326  output_field_ptr->set_fe_data(feo->dh(), feo->mapping<1>(), feo->mapping<2>(), feo->mapping<3>(), &output_vec[sbi]);
327  data_.output_field[sbi].set_field(Model::mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
328  }
329 
330  // set time marks for writing the output
331  data_.output_fields.initialize(Model::output_stream_, input_rec.val<Input::Record>("output"), this->time());
332 
333  // allocate matrix and vector structures
334  ls = new LinSys*[Model::n_substances()];
335  ls_dt = new LinSys*[Model::n_substances()];
336  solution_elem_ = new double*[Model::n_substances()];
337  for (unsigned int sbi = 0; sbi < Model::n_substances(); sbi++) {
338  ls[sbi] = new LinSys_PETSC(feo->dh()->distr());
339  ( (LinSys_PETSC *)ls[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
340  ls[sbi]->set_solution(NULL);
341 
342  ls_dt[sbi] = new LinSys_PETSC(feo->dh()->distr());
343  ( (LinSys_PETSC *)ls_dt[sbi] )->set_from_input( input_rec.val<Input::Record>("solver") );
344  solution_elem_[sbi] = new double[Model::mesh_->get_el_ds()->lsize()];
345  }
346  stiffness_matrix = new Mat[Model::n_substances()];
347  mass_matrix = new Mat[Model::n_substances()];
348  rhs = new Vec[Model::n_substances()];
349  mass_vec = new Vec[Model::n_substances()];
350 
351 
352  // initialization of balance object
353  if (Model::balance_ != nullptr)
354  {
355  Model::balance_->allocate(feo->dh()->distr()->lsize(),
356  max(feo->fe<1>()->n_dofs(), max(feo->fe<2>()->n_dofs(), feo->fe<3>()->n_dofs())));
357  }
358 
359 }
360 
361 
362 template<class Model>
364 {
365  delete Model::time_;
366 
367  if (gamma.size() > 0) {
368  // initialize called
369 
370  for (auto &vec : output_vec) VecDestroy(&vec);
371 
372  for (unsigned int i=0; i<Model::n_substances(); i++)
373  {
374  delete ls[i];
375  delete[] solution_elem_[i];
376  delete ls_dt[i];
377  MatDestroy(&stiffness_matrix[i]);
378  MatDestroy(&mass_matrix[i]);
379  VecDestroy(&rhs[i]);
380  VecDestroy(&mass_vec[i]);
381  }
382  delete[] ls;
383  delete[] solution_elem_;
384  delete[] ls_dt;
385  delete[] stiffness_matrix;
386  delete[] mass_matrix;
387  delete[] rhs;
388  delete[] mass_vec;
389  delete feo;
390 
391  }
392 
393 }
394 
395 
396 template<class Model>
398 {
399  VecScatter output_scatter;
400  VecScatterCreateToZero(ls[0]->get_solution(), &output_scatter, PETSC_NULL);
401  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
402  {
403  // gather solution to output_vec[sbi]
404  VecScatterBegin(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
405  VecScatterEnd(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
406  }
407  VecScatterDestroy(&(output_scatter));
408 
409 }
410 
411 
412 
413 template<class Model>
415 {
416  START_TIMER(Model::ModelEqData::name());
417  data_.mark_input_times( *(Model::time_) );
418  data_.set_time(Model::time_->step(), LimitSide::left);
419 
420 
421  // set initial conditions
423  for (unsigned int sbi = 0; sbi < Model::n_substances(); sbi++)
424  ( (LinSys_PETSC *)ls[sbi] )->set_initial_guess_nonzero();
425 
426  // check first time assembly - needs preallocation
428 
429  // after preallocation we assemble the matrices and vectors required for mass balance
430  if (Model::balance_ != nullptr)
431  {
432  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
433  {
434  Model::balance_->calculate_mass(Model::subst_idx[sbi], ls[sbi]->get_solution());
435  Model::balance_->calculate_source(Model::subst_idx[sbi], ls[sbi]->get_solution());
436  Model::balance_->calculate_flux(Model::subst_idx[sbi], ls[sbi]->get_solution());
437 
438  // add sources due to sorption
439  vector<double> masses(Model::mesh_->region_db().bulk_size());
440  double mass = 0;
441  Model::balance_->calculate_mass(Model::subst_idx[sbi], ls[sbi]->get_solution(), masses);
442  for (auto reg_mass : masses)
443  mass += reg_mass;
444 
445  sorption_sources[sbi] = mass;
446  }
447  }
448 
449  output_data();
450 }
451 
452 
453 template<class Model>
455 {
456  // preallocate system matrix
457  for (unsigned int i=0; i<Model::n_substances(); i++)
458  {
459  // preallocate system matrix
460  ls[i]->start_allocation();
461  stiffness_matrix[i] = NULL;
462  rhs[i] = NULL;
463 
464  // preallocate mass matrix
465  ls_dt[i]->start_allocation();
466  mass_matrix[i] = NULL;
467  }
470  set_sources();
472 
473  allocation_done = true;
474 }
475 
476 
477 
478 template<class Model>
480 {
481  START_TIMER("DG-ONE STEP");
482 
483  Model::time_->next_time();
484  Model::time_->view("TDG");
485 
486  START_TIMER("data reinit");
487  data_.set_time(Model::time_->step(), LimitSide::left);
488  END_TIMER("data reinit");
489 
490  // assemble mass matrix
491  if (mass_matrix[0] == NULL || data_.subset(FieldFlag::in_time_term).changed() )
492  {
493  for (unsigned int i=0; i<Model::n_substances(); i++)
494  {
496  ls_dt[i]->mat_zero_entries();
497  }
499  for (unsigned int i=0; i<Model::n_substances(); i++)
500  {
501  ls_dt[i]->finish_assembly();
502  // construct mass_vec for initial time
503  if (mass_matrix[i] == NULL)
504  {
505  VecDuplicate(ls[i]->get_solution(), &mass_vec[i]);
506  MatMult(*(ls_dt[i]->get_matrix()), ls[i]->get_solution(), mass_vec[i]);
507  if (Model::balance_ != nullptr && Model::balance_->cumulative())
508  {
509  double total_mass = 0;
510  VecSum(mass_vec[i], &total_mass);
511  sorption_sources[i] -= total_mass;
512  }
513  MatConvert(*( ls_dt[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &mass_matrix[i]);
514  }
515  else
516  MatCopy(*( ls_dt[i]->get_matrix() ), mass_matrix[i], DIFFERENT_NONZERO_PATTERN);
517  }
518  }
519 
520  // assemble stiffness matrix
521  if (stiffness_matrix[0] == NULL
522  || data_.subset(FieldFlag::in_main_matrix).changed()
523  || Model::flux_changed)
524  {
525  // new fluxes can change the location of Neumann boundary,
526  // thus stiffness matrix must be reassembled
527  for (unsigned int i=0; i<Model::n_substances(); i++)
528  {
529  ls[i]->start_add_assembly();
530  ls[i]->mat_zero_entries();
531  }
533  for (unsigned int i=0; i<Model::n_substances(); i++)
534  {
535  ls[i]->finish_assembly();
536 
537  if (stiffness_matrix[i] == NULL)
538  MatConvert(*( ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
539  else
540  MatCopy(*( ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
541  }
542  }
543 
544  // assemble right hand side (due to sources and boundary conditions)
545  if (rhs[0] == NULL
546  || data_.subset(FieldFlag::in_rhs).changed()
547  || Model::flux_changed)
548  {
549  for (unsigned int i=0; i<Model::n_substances(); i++)
550  {
551  ls[i]->start_add_assembly();
552  ls[i]->rhs_zero_entries();
553  }
554  set_sources();
556  for (unsigned int i=0; i<Model::n_substances(); i++)
557  {
558  ls[i]->finish_assembly();
559 
560  if (rhs[i] == nullptr) VecDuplicate(*( ls[i]->get_rhs() ), &rhs[i]);
561  VecCopy(*( ls[i]->get_rhs() ), rhs[i]);
562  }
563  }
564 
565  Model::flux_changed = false;
566 
567 
568  /* Apply backward Euler time integration.
569  *
570  * Denoting A the stiffness matrix and M the mass matrix, the algebraic system at the k-th time level reads
571  *
572  * (1/dt M + A)u^k = f + 1/dt M.u^{k-1}
573  *
574  * Hence we modify at each time level the right hand side:
575  *
576  * f^k = f + 1/dt M u^{k-1},
577  *
578  * where f stands for the term stemming from the force and boundary conditions.
579  * Accordingly, we set
580  *
581  * A^k = A + 1/dt M.
582  *
583  */
584  Mat m;
585  START_TIMER("solve");
586  for (unsigned int i=0; i<Model::n_substances(); i++)
587  {
588  MatConvert(stiffness_matrix[i], MATSAME, MAT_INITIAL_MATRIX, &m);
589  MatAXPY(m, 1./Model::time_->dt(), mass_matrix[i], SUBSET_NONZERO_PATTERN);
590  ls[i]->set_matrix(m, DIFFERENT_NONZERO_PATTERN);
591  Vec w;
592  VecDuplicate(rhs[i], &w);
593  VecWAXPY(w, 1./Model::time_->dt(), mass_vec[i], rhs[i]);
594  ls[i]->set_rhs(w);
595 
596  VecDestroy(&w);
597  MatDestroy(&m);
598 
599  ls[i]->solve();
600 
601  // update mass_vec due to possible changes in mass matrix
602  MatMult(*(ls_dt[i]->get_matrix()), ls[i]->get_solution(), mass_vec[i]);
603  }
604  END_TIMER("solve");
605 
607 
608  END_TIMER("DG-ONE STEP");
609 }
610 
611 
612 template<class Model>
614 {
615  // calculate element averages of solution
616  for (unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
617  {
618  typename DOFHandlerBase::CellIterator elem = Model::mesh_->element(feo->dh()->el_index(i_cell));
619 
620  unsigned int n_dofs;
621  switch (elem->dim())
622  {
623  case 1:
624  n_dofs = feo->fe<1>()->n_dofs();
625  break;
626  case 2:
627  n_dofs = feo->fe<2>()->n_dofs();
628  break;
629  case 3:
630  n_dofs = feo->fe<3>()->n_dofs();
631  break;
632  }
633 
634  unsigned int dof_indices[n_dofs];
635  feo->dh()->get_dof_indices(elem, dof_indices);
636 
637  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
638  {
639  solution_elem_[sbi][i_cell] = 0;
640 
641  for (unsigned int j=0; j<n_dofs; ++j)
642  solution_elem_[sbi][i_cell] += ls[sbi]->get_solution_array()[dof_indices[j]-feo->dh()->distr()->begin()];
643 
644  solution_elem_[sbi][i_cell] = max(solution_elem_[sbi][i_cell]/n_dofs, 0.);
645  }
646  }
647 }
648 
649 
650 
651 
652 template<class Model>
654 {
655  //if (!Model::time_->is_current( Model::time_->marks().type_output() )) return;
656 
657 
658  START_TIMER("DG-OUTPUT");
659 
660  // gather the solution from all processors
661  data_.output_fields.set_time( this->time().step(), LimitSide::left);
662  //if (data_.output_fields.is_field_output_time(data_.output_field, this->time().step()) )
664  data_.output_fields.output(this->time().step());
665 
666  Model::output_data();
667 
668  END_TIMER("DG-OUTPUT");
669 }
670 
671 
672 template<class Model>
674 {
675  if (Model::balance_ != nullptr && Model::balance_->cumulative())
676  {
677  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
678  {
679  Model::balance_->calculate_cumulative_sources(Model::subst_idx[sbi], ls[sbi]->get_solution(), Model::time_->dt());
680  Model::balance_->calculate_cumulative_fluxes(Model::subst_idx[sbi], ls[sbi]->get_solution(), Model::time_->dt());
681 
682  // add sources due to sorption
683  vector<double> masses(Model::mesh_->region_db().bulk_size());
684  double mass = 0;
685  Model::balance_->calculate_mass(Model::subst_idx[sbi], ls[sbi]->get_solution(), masses);
686  for (auto reg_mass : masses)
687  mass += reg_mass;
688  double total_mass = 0;
689  VecSum(mass_vec[sbi], &total_mass);
690 
691  Model::balance_->add_cumulative_source(Model::subst_idx[sbi], mass-total_mass-sorption_sources[sbi]);
692  sorption_sources[sbi] = mass - total_mass;
693  }
694  }
695 }
696 
697 
698 template<class Model>
700 {
701  if (Model::balance_ != nullptr)
702  {
703  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
704  {
705  Model::balance_->calculate_mass(Model::subst_idx[sbi], ls[sbi]->get_solution());
706  Model::balance_->calculate_source(Model::subst_idx[sbi], ls[sbi]->get_solution());
707  Model::balance_->calculate_flux(Model::subst_idx[sbi], ls[sbi]->get_solution());
708  }
709  }
710 }
711 
712 
713 template<class Model>
715 {
716  START_TIMER("assemble_mass");
717  if (Model::balance_ != nullptr)
718  Model::balance_->start_mass_assembly(Model::subst_idx);
719  assemble_mass_matrix<1>();
720  assemble_mass_matrix<2>();
721  assemble_mass_matrix<3>();
722  if (Model::balance_ != nullptr)
723  Model::balance_->finish_mass_assembly(Model::subst_idx);
724  END_TIMER("assemble_mass");
725 }
726 
727 
728 template<class Model> template<unsigned int dim>
730 {
731  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(), update_values | update_JxW_values | update_quadrature_points);
732  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
733  vector<int> dof_indices(ndofs);
734  PetscScalar local_mass_matrix[ndofs*ndofs];
735  vector<PetscScalar> local_mass_balance_vector(ndofs);
736 
737  // assemble integral over elements
738  for (unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
739  {
740  typename DOFHandlerBase::CellIterator cell = Model::mesh_->element(feo->dh()->el_index(i_cell));
741  if (cell->dim() != dim) continue;
742 
743  fe_values.reinit(cell);
744  feo->dh()->get_dof_indices(cell, (unsigned int*)&(dof_indices[0]));
745  ElementAccessor<3> ele_acc = cell->element_accessor();
746 
747  Model::compute_mass_matrix_coefficient(fe_values.point_list(), ele_acc, mm_coef);
748  Model::compute_retardation_coefficient(fe_values.point_list(), ele_acc, ret_coef);
749 
750  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
751  {
752  // assemble the local mass matrix
753  for (unsigned int i=0; i<ndofs; i++)
754  {
755  for (unsigned int j=0; j<ndofs; j++)
756  {
757  local_mass_matrix[i*ndofs+j] = 0;
758  for (unsigned int k=0; k<qsize; k++)
759  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);
760  }
761  }
762 
763  if (Model::balance_ != nullptr)
764  {
765  for (unsigned int i=0; i<ndofs; i++)
766  {
767  local_mass_balance_vector[i] = 0;
768  for (unsigned int k=0; k<qsize; k++)
769  local_mass_balance_vector[i] += mm_coef[k]*fe_values.shape_value(i,k)*fe_values.JxW(k);
770  }
771  }
772 
773  Model::balance_->add_mass_matrix_values(Model::subst_idx[sbi], ele_acc.region().bulk_idx(), dof_indices, local_mass_balance_vector);
774  ls_dt[sbi]->mat_set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), local_mass_matrix);
775  }
776  }
777 }
778 
779 
780 
781 
782 template<class Model>
784 {
785  START_TIMER("assemble_stiffness");
786  START_TIMER("assemble_volume_integrals");
787  assemble_volume_integrals<1>();
788  assemble_volume_integrals<2>();
789  assemble_volume_integrals<3>();
790  END_TIMER("assemble_volume_integrals");
791 
792  START_TIMER("assemble_fluxes_boundary");
793  assemble_fluxes_boundary<1>();
794  assemble_fluxes_boundary<2>();
795  assemble_fluxes_boundary<3>();
796  END_TIMER("assemble_fluxes_boundary");
797 
798  START_TIMER("assemble_fluxes_elem_elem");
799  assemble_fluxes_element_element<1>();
800  assemble_fluxes_element_element<2>();
801  assemble_fluxes_element_element<3>();
802  END_TIMER("assemble_fluxes_elem_elem");
803 
804  START_TIMER("assemble_fluxes_elem_side");
805  assemble_fluxes_element_side<1>();
806  assemble_fluxes_element_side<2>();
807  assemble_fluxes_element_side<3>();
808  END_TIMER("assemble_fluxes_elem_side");
809  END_TIMER("assemble_stiffness");
810 }
811 
812 
813 
814 template<class Model>
815 template<unsigned int dim>
817 {
818  FEValues<dim,3> fv_rt(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe_rt<dim>(),
820  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
822  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
823  unsigned int dof_indices[ndofs];
824  vector<arma::vec3> velocity(qsize);
825  vector<arma::vec> sources_sigma(qsize, arma::vec(Model::n_substances()));
826  PetscScalar local_matrix[ndofs*ndofs];
827 
828  // assemble integral over elements
829  for (unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
830  {
831  typename DOFHandlerBase::CellIterator cell = Model::mesh_->element(feo->dh()->el_index(i_cell));
832  if (cell->dim() != dim) continue;
833 
834  fe_values.reinit(cell);
835  fv_rt.reinit(cell);
836  ElementAccessor<3> ele_acc = cell->element_accessor();
837  feo->dh()->get_dof_indices(cell, dof_indices);
838 
839  calculate_velocity(cell, velocity, fv_rt);
840  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), velocity, ele_acc, ad_coef, dif_coef);
841  Model::compute_sources_sigma(fe_values.point_list(), ele_acc, sources_sigma);
842 
843  // assemble the local stiffness matrix
844  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
845  {
846  for (unsigned int i=0; i<ndofs; i++)
847  for (unsigned int j=0; j<ndofs; j++)
848  local_matrix[i*ndofs+j] = 0;
849 
850  for (unsigned int k=0; k<qsize; k++)
851  {
852  for (unsigned int i=0; i<ndofs; i++)
853  {
854  arma::vec3 Kt_grad_i = dif_coef[sbi][k].t()*fe_values.shape_grad(i,k);
855  double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.shape_grad(i,k));
856 
857  for (unsigned int j=0; j<ndofs; j++)
858  local_matrix[i*ndofs+j] += (arma::dot(Kt_grad_i, fe_values.shape_grad(j,k))
859  -fe_values.shape_value(j,k)*ad_dot_grad_i
860  +sources_sigma[k][sbi]*fe_values.shape_value(j,k)*fe_values.shape_value(i,k))*fe_values.JxW(k);
861  }
862  }
863  ls[sbi]->mat_set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, local_matrix);
864  }
865  }
866 }
867 
868 
869 template<class Model>
871 {
872  START_TIMER("assemble_sources");
873  if (Model::balance_ != nullptr)
874  Model::balance_->start_source_assembly(Model::subst_idx);
875  set_sources<1>();
876  set_sources<2>();
877  set_sources<3>();
878  if (Model::balance_ != nullptr)
879  Model::balance_->finish_source_assembly(Model::subst_idx);
880  END_TIMER("assemble_sources");
881 }
882 
883 template<class Model>
884 template<unsigned int dim>
886 {
887  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
889  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
890  vector<arma::vec> sources_conc(qsize, arma::vec(Model::n_substances())),
891  sources_density(qsize, arma::vec(Model::n_substances())),
892  sources_sigma(qsize, arma::vec(Model::n_substances()));
893  vector<int> dof_indices(ndofs);
894  PetscScalar local_rhs[ndofs];
895  vector<PetscScalar> local_source_balance_vector(ndofs), local_source_balance_rhs(ndofs);
896  double source;
897 
898  // assemble integral over elements
899  for (unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
900  {
901  typename DOFHandlerBase::CellIterator cell = Model::mesh_->element(feo->dh()->el_index(i_cell));
902  if (cell->dim() != dim) continue;
903 
904  fe_values.reinit(cell);
905  feo->dh()->get_dof_indices(cell, (unsigned int *)&(dof_indices[0]));
906 
907  Model::compute_source_coefficients(fe_values.point_list(), cell->element_accessor(), sources_conc, sources_density, sources_sigma);
908 
909  // assemble the local stiffness matrix
910  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
911  {
912  fill_n(local_rhs, ndofs, 0);
913  local_source_balance_vector.assign(ndofs, 0);
914  local_source_balance_rhs.assign(ndofs, 0);
915 
916  // compute sources
917  for (unsigned int k=0; k<qsize; k++)
918  {
919  source = (sources_density[k][sbi] + sources_conc[k][sbi]*sources_sigma[k][sbi])*fe_values.JxW(k);
920 
921  for (unsigned int i=0; i<ndofs; i++)
922  local_rhs[i] += source*fe_values.shape_value(i,k);
923  }
924  ls[sbi]->rhs_set_values(ndofs, &(dof_indices[0]), local_rhs);
925 
926  if (Model::balance_ != nullptr)
927  {
928  for (unsigned int i=0; i<ndofs; i++)
929  {
930  for (unsigned int k=0; k<qsize; k++)
931  local_source_balance_vector[i] -= sources_sigma[k][sbi]*fe_values.shape_value(i,k)*fe_values.JxW(k);
932 
933  local_source_balance_rhs[i] += local_rhs[i];
934  }
935  Model::balance_->add_source_matrix_values(Model::subst_idx[sbi], cell->region().bulk_idx(), dof_indices, local_source_balance_vector);
936  Model::balance_->add_source_vec_values(Model::subst_idx[sbi], cell->region().bulk_idx(), dof_indices, local_source_balance_rhs);
937  }
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<arma::vec> 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++)
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] = data_.dg_penalty.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  arma::vec dg_penalty;
1090  double gamma_l;
1091 
1092  // assemble boundary integral
1093  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1094  {
1095  Edge *edg = &Model::mesh_->edges[feo->dh()->edge_index(iedg)];
1096  if (edg->n_sides > 1) continue;
1097  // check spatial dimension
1098  if (edg->side(0)->dim() != dim-1) continue;
1099  // skip edges lying not on the boundary
1100  if (edg->side(0)->cond() == NULL) continue;
1101 
1102  SideIter side = edg->side(0);
1103  typename DOFHandlerBase::CellIterator cell = side->element();
1104  ElementAccessor<3> ele_acc = cell->element_accessor();
1105  feo->dh()->get_dof_indices(cell, side_dof_indices);
1106  fe_values_side.reinit(cell, side->el_idx());
1107  fsv_rt.reinit(cell, side->el_idx());
1108 
1109  calculate_velocity(cell, side_velocity, fsv_rt);
1110  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), side_velocity, ele_acc, ad_coef, dif_coef);
1111  dg_penalty = data_.dg_penalty.value(cell->centre(), ele_acc);
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), dg_penalty[sbi], 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<arma::vec> frac_sigma(qsize, arma::vec(Model::n_substances()));
1199  vector<double> csection_lower(qsize), csection_higher(qsize), por_lower(qsize), por_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  data_.fracture_sigma.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), frac_sigma);
1239  data_.porosity.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), por_lower);
1240  data_.porosity.value_list(fe_values_vb.point_list(), cell->element_accessor(), por_higher);
1241 
1242  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1243  {
1244  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
1245  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
1246  local_matrix[i*(n_dofs[0]+n_dofs[1])+j] = 0;
1247 
1248  // set transmission conditions
1249  for (unsigned int k=0; k<qsize; k++)
1250  {
1251  /* The communication flux has two parts:
1252  * - "diffusive" term containing sigma
1253  * - "advective" term representing usual upwind
1254  *
1255  * The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
1256  * than b and A in the manual.
1257  * In calculation of sigma there appears one more csection_lower in the denominator.
1258  */
1259  double sigma = frac_sigma[k][sbi]*arma::dot(dif_coef_edg[0][sbi][k]*fe_values_side.normal_vector(k),fe_values_side.normal_vector(k))*
1260  2*csection_higher[k]*csection_higher[k]/(csection_lower[k]*csection_lower[k]);
1261 
1262  double transport_flux = arma::dot(ad_coef_edg[1][sbi][k], fe_values_side.normal_vector(k));
1263  double por_lower_over_higher = por_lower[k]/por_higher[k];
1264 
1265  comm_flux[0][0] = (sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1266  comm_flux[0][1] = -(sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1267  comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1268  comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1269 
1270  for (int n=0; n<2; n++)
1271  {
1272  if (!feo->dh()->el_is_local(element_id[n])) continue;
1273 
1274  for (unsigned int i=0; i<n_dofs[n]; i++)
1275  for (int m=0; m<2; m++)
1276  for (unsigned int j=0; j<n_dofs[m]; j++)
1277  local_matrix[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
1278  comm_flux[m][n]*fv_sb[m]->shape_value(j,k)*fv_sb[n]->shape_value(i,k);
1279  }
1280  }
1281  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);
1282  }
1283  }
1284 
1285 }
1286 
1287 
1288 
1289 
1290 
1291 
1292 template<class Model>
1294 {
1295  START_TIMER("assemble_bc");
1296  if (Model::balance_ != nullptr)
1297  Model::balance_->start_flux_assembly(Model::subst_idx);
1298  set_boundary_conditions<1>();
1299  set_boundary_conditions<2>();
1300  set_boundary_conditions<3>();
1301  if (Model::balance_ != nullptr)
1302  Model::balance_->finish_flux_assembly(Model::subst_idx);
1303  END_TIMER("assemble_bc");
1304 }
1305 
1306 
1307 template<class Model>
1308 template<unsigned int dim>
1310 {
1311  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1313  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1314  update_values);
1315  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1316  vector<int> side_dof_indices(ndofs);
1317  unsigned int loc_b=0;
1318  double local_rhs[ndofs];
1319  vector<PetscScalar> local_flux_balance_vector(ndofs);
1320  PetscScalar local_flux_balance_rhs;
1321  vector<arma::vec> bc_values(qsize, arma::vec(Model::n_substances()));
1322  vector<double> bc_fluxes(qsize),
1323  bc_sigma(qsize),
1324  bc_ref_values(qsize);
1325  vector<double> csection(qsize);
1326  vector<arma::vec3> velocity;
1327 
1328  for (unsigned int loc_el = 0; loc_el < Model::mesh_->get_el_ds()->lsize(); loc_el++)
1329  {
1330  ElementFullIter elm = Model::mesh_->element(feo->dh()->el_index(loc_el));
1331  if (elm->boundary_idx_ == nullptr) continue;
1332 
1333  FOR_ELEMENT_SIDES(elm,si)
1334  {
1335  Edge *edg = elm->side(si)->edge();
1336  if (edg->n_sides > 1) continue;
1337  // skip edges lying not on the boundary
1338  if (edg->side(0)->cond() == NULL) continue;
1339 
1340  if (edg->side(0)->dim() != dim-1)
1341  {
1342  if (edg->side(0)->cond() != nullptr) ++loc_b;
1343  continue;
1344  }
1345 
1346  SideIter side = edg->side(0);
1347  typename DOFHandlerBase::CellIterator cell = Model::mesh_->element.full_iter(side->element());
1348  ElementAccessor<3> ele_acc = side->cond()->element_accessor();
1349 
1350  arma::uvec bc_type;
1351  Model::get_bc_type(ele_acc, bc_type);
1352 
1353  fe_values_side.reinit(cell, side->el_idx());
1354  fsv_rt.reinit(cell, side->el_idx());
1355  calculate_velocity(cell, velocity, fsv_rt);
1356 
1357  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), velocity, side->element()->element_accessor(), ad_coef, dif_coef);
1358  data_.cross_section.value_list(fe_values_side.point_list(), side->element()->element_accessor(), csection);
1359  // The b.c. data are fetched for all possible b.c. types since we allow
1360  // different bc_type for each substance.
1361  data_.bc_dirichlet_value.value_list(fe_values_side.point_list(), ele_acc, bc_values);
1362 
1363  feo->dh()->get_dof_indices(cell, (unsigned int *)&(side_dof_indices[0]));
1364 
1365  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1366  {
1367  fill_n(local_rhs, ndofs, 0);
1368  local_flux_balance_vector.assign(ndofs, 0);
1369  local_flux_balance_rhs = 0;
1370 
1371  double side_flux = 0;
1372  for (unsigned int k=0; k<qsize; k++)
1373  side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
1374  double transport_flux = side_flux/side->measure();
1375 
1376  if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
1377  {
1378  for (unsigned int k=0; k<qsize; k++)
1379  {
1380  double bc_term = -transport_flux*bc_values[k][sbi]*fe_values_side.JxW(k);
1381  for (unsigned int i=0; i<ndofs; i++)
1382  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k);
1383  }
1384  if (Model::balance_ != nullptr)
1385  for (unsigned int i=0; i<ndofs; i++)
1386  local_flux_balance_rhs -= local_rhs[i];
1387  }
1388  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_dirichlet)
1389  {
1390  for (unsigned int k=0; k<qsize; k++)
1391  {
1392  double bc_term = gamma[sbi][side->cond_idx()]*bc_values[k][sbi]*fe_values_side.JxW(k);
1393  arma::vec3 bc_grad = -bc_values[k][sbi]*fe_values_side.JxW(k)*dg_variant*(arma::trans(dif_coef[sbi][k])*fe_values_side.normal_vector(k));
1394  for (unsigned int i=0; i<ndofs; i++)
1395  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k)
1396  + arma::dot(bc_grad,fe_values_side.shape_grad(i,k));
1397  }
1398  if (Model::balance_ != nullptr)
1399  {
1400  for (unsigned int k=0; k<qsize; k++)
1401  {
1402  for (unsigned int i=0; i<ndofs; i++)
1403  {
1404  local_flux_balance_vector[i] += (arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.shape_value(i,k)
1405  - arma::dot(dif_coef[sbi][k]*fe_values_side.shape_grad(i,k),fe_values_side.normal_vector(k))
1406  + gamma[sbi][side->cond_idx()]*fe_values_side.shape_value(i,k))*fe_values_side.JxW(k);
1407  }
1408  }
1409  if (Model::time_->tlevel() > 0)
1410  for (unsigned int i=0; i<ndofs; i++)
1411  local_flux_balance_rhs -= local_rhs[i];
1412  }
1413  }
1414  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_total_flux)
1415  {
1416  Model::get_flux_bc_data(sbi, fe_values_side.point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1417  for (unsigned int k=0; k<qsize; k++)
1418  {
1419  double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.JxW(k);
1420  for (unsigned int i=0; i<ndofs; i++)
1421  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k);
1422  }
1423 
1424  if (Model::balance_ != nullptr)
1425  {
1426  for (unsigned int i=0; i<ndofs; i++)
1427  {
1428  for (unsigned int k=0; k<qsize; k++)
1429  local_flux_balance_vector[i] += csection[k]*bc_sigma[k]*fe_values_side.JxW(k)*fe_values_side.shape_value(i,k);
1430  local_flux_balance_rhs -= local_rhs[i];
1431  }
1432  }
1433  }
1434  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_diffusive_flux)
1435  {
1436  Model::get_flux_bc_data(sbi, fe_values_side.point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1437  for (unsigned int k=0; k<qsize; k++)
1438  {
1439  double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.JxW(k);
1440  for (unsigned int i=0; i<ndofs; i++)
1441  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k);
1442  }
1443 
1444  if (Model::balance_ != nullptr)
1445  {
1446  for (unsigned int i=0; i<ndofs; i++)
1447  {
1448  for (unsigned int k=0; k<qsize; k++)
1449  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);
1450  local_flux_balance_rhs -= local_rhs[i];
1451  }
1452  }
1453  }
1454  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux >= 0)
1455  {
1456  if (Model::balance_ != nullptr)
1457  {
1458  for (unsigned int k=0; k<qsize; k++)
1459  {
1460  for (unsigned int i=0; i<ndofs; i++)
1461  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);
1462  }
1463  }
1464  }
1465  ls[sbi]->rhs_set_values(ndofs, &(side_dof_indices[0]), local_rhs);
1466 
1467  if (Model::balance_ != nullptr)
1468  {
1469  Model::balance_->add_flux_matrix_values(Model::subst_idx[sbi], loc_b, side_dof_indices, local_flux_balance_vector);
1470  Model::balance_->add_flux_vec_value(Model::subst_idx[sbi], loc_b, local_flux_balance_rhs);
1471  }
1472  }
1473  ++loc_b;
1474  }
1475  }
1476 }
1477 
1478 
1479 
1480 template<class Model>
1481 template<unsigned int dim>
1483 {
1484  OLD_ASSERT(cell->dim() == dim, "Element dimension mismatch!");
1485 
1486  velocity.resize(fv.n_points());
1487 
1488  for (unsigned int k=0; k<fv.n_points(); k++)
1489  {
1490  velocity[k].zeros();
1491  for (unsigned int sid=0; sid<cell->n_sides(); sid++)
1492  velocity[k] += fv.shape_vector(sid,k) * Model::mh_dh->side_flux( *(cell->side(sid)) );
1493  }
1494 }
1495 
1496 
1497 
1498 
1499 
1500 template<class Model>
1502  const int s1,
1503  const int s2,
1504  const int K_size,
1505  const vector<arma::mat33> &K1,
1506  const vector<arma::mat33> &K2,
1507  const vector<double> &fluxes,
1508  const arma::vec3 &normal_vector,
1509  const double alpha1,
1510  const double alpha2,
1511  double &gamma,
1512  double *omega,
1513  double &transport_flux)
1514 {
1515  double delta[2];
1516  double h = 0;
1517  double local_alpha = 0;
1518 
1519  OLD_ASSERT(edg.side(s1)->valid(), "Invalid side of an edge.");
1520  SideIter s = edg.side(s1);
1521 
1522  // calculate the side diameter
1523  if (s->dim() == 0)
1524  {
1525  h = 1;
1526  }
1527  else
1528  {
1529  for (unsigned int i=0; i<s->n_nodes(); i++)
1530  for (unsigned int j=i+1; j<s->n_nodes(); j++)
1531  h = max(h, s->node(i)->distance(*s->node(j)));
1532  }
1533 
1534  // calculate the total in- and out-flux through the edge
1535  double pflux = 0, nflux = 0;
1536  for (int i=0; i<edg.n_sides; i++)
1537  {
1538  if (fluxes[i] > 0)
1539  pflux += fluxes[i];
1540  else
1541  nflux += fluxes[i];
1542  }
1543 
1544  // calculate the flux from s1 to s2
1545  if (fluxes[s2] > 0 && fluxes[s1] < 0 && s1 < s2)
1546  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
1547  else if (fluxes[s2] < 0 && fluxes[s1] > 0 && s1 < s2)
1548  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
1549  else if (s1==s2)
1550  transport_flux = fluxes[s1];
1551  else
1552  transport_flux = 0;
1553 
1554  gamma = 0.5*fabs(transport_flux);
1555 
1556 
1557  // determine local DG penalty
1558  local_alpha = max(alpha1, alpha2);
1559 
1560  if (s1 == s2)
1561  {
1562  omega[0] = 1;
1563 
1564  // delta is set to the average value of Kn.n on the side
1565  delta[0] = 0;
1566  for (int k=0; k<K_size; k++)
1567  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1568  delta[0] /= K_size;
1569 
1570  gamma += local_alpha/h*(delta[0]);
1571  }
1572  else
1573  {
1574  delta[0] = 0;
1575  delta[1] = 0;
1576  for (int k=0; k<K_size; k++)
1577  {
1578  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1579  delta[1] += dot(K2[k]*normal_vector,normal_vector);
1580  }
1581  delta[0] /= K_size;
1582  delta[1] /= K_size;
1583 
1584  double delta_sum = delta[0] + delta[1];
1585 
1586  if (delta_sum > numeric_limits<double>::epsilon())
1587  {
1588  omega[0] = delta[1]/delta_sum;
1589  omega[1] = delta[0]/delta_sum;
1590  gamma += local_alpha/h*(delta[0]*delta[1]/delta_sum);
1591  }
1592  else
1593  for (int i=0; i<2; i++) omega[i] = 0;
1594  }
1595 }
1596 
1597 
1598 
1599 
1600 
1601 
1602 template<class Model>
1604  const int K_size,
1605  const vector<arma::mat33> &K,
1606  const double flux,
1607  const arma::vec3 &normal_vector,
1608  const double alpha,
1609  double &gamma)
1610 {
1611  double delta = 0, h = 0;
1612 
1613  // calculate the side diameter
1614  if (side->dim() == 0)
1615  {
1616  h = 1;
1617  }
1618  else
1619  {
1620  for (unsigned int i=0; i<side->n_nodes(); i++)
1621  for (unsigned int j=i+1; j<side->n_nodes(); j++)
1622  h = max(h, side->node(i)->distance( *side->node(j) ));
1623  }
1624 
1625  // delta is set to the average value of Kn.n on the side
1626  for (int k=0; k<K_size; k++)
1627  delta += dot(K[k]*normal_vector,normal_vector);
1628  delta /= K_size;
1629 
1630  gamma = 0.5*fabs(flux) + alpha/h*delta;
1631 }
1632 
1633 
1634 
1635 
1636 
1637 template<class Model>
1639 {
1640  START_TIMER("set_init_cond");
1641  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1642  ls[sbi]->start_allocation();
1643  prepare_initial_condition<1>();
1644  prepare_initial_condition<2>();
1645  prepare_initial_condition<3>();
1646 
1647  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1648  ls[sbi]->start_add_assembly();
1649  prepare_initial_condition<1>();
1650  prepare_initial_condition<2>();
1651  prepare_initial_condition<3>();
1652 
1653  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1654  {
1655  ls[sbi]->finish_assembly();
1656  ls[sbi]->solve();
1657  }
1658  END_TIMER("set_init_cond");
1659 }
1660 
1661 template<class Model>
1662 template<unsigned int dim>
1664 {
1665  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1667  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1668  unsigned int dof_indices[ndofs];
1669  double matrix[ndofs*ndofs], rhs[ndofs];
1670  std::vector<arma::vec> init_values(qsize);
1671 
1672  for (unsigned int k=0; k<qsize; k++)
1673  init_values[k].resize(Model::n_substances());
1674 
1675  for (unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
1676  {
1677  typename DOFHandlerBase::CellIterator elem = Model::mesh_->element(feo->dh()->el_index(i_cell));
1678  if (elem->dim() != dim) continue;
1679 
1680  ElementAccessor<3> ele_acc = elem->element_accessor();
1681  feo->dh()->get_dof_indices(elem, dof_indices);
1682  fe_values.reinit(elem);
1683 
1684  Model::compute_init_cond(fe_values.point_list(), ele_acc, init_values);
1685 
1686  for (unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1687  {
1688  for (unsigned int i=0; i<ndofs; i++)
1689  {
1690  rhs[i] = 0;
1691  for (unsigned int j=0; j<ndofs; j++)
1692  matrix[i*ndofs+j] = 0;
1693  }
1694 
1695  for (unsigned int k=0; k<qsize; k++)
1696  {
1697  double rhs_term = init_values[k](sbi)*fe_values.JxW(k);
1698 
1699  for (unsigned int i=0; i<ndofs; i++)
1700  {
1701  for (unsigned int j=0; j<ndofs; j++)
1702  matrix[i*ndofs+j] += fe_values.shape_value(i,k)*fe_values.shape_value(j,k)*fe_values.JxW(k);
1703 
1704  rhs[i] += fe_values.shape_value(i,k)*rhs_term;
1705  }
1706  }
1707  ls[sbi]->set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, matrix, rhs);
1708  }
1709  }
1710 }
1711 
1712 
1713 template<class Model>
1714 void TransportDG<Model>::get_par_info(int * &el_4_loc, Distribution * &el_ds)
1715 {
1716  el_4_loc = Model::mesh_->get_el_4_loc();
1717  el_ds = Model::mesh_->get_el_ds();
1718 }
1719 
1720 
1721 template<class Model>
1723 {
1724  if (solution_changed)
1725  {
1726  for (unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
1727  {
1728  typename DOFHandlerBase::CellIterator elem = Model::mesh_->element(feo->dh()->el_index(i_cell));
1729 
1730  unsigned int n_dofs;
1731  switch (elem->dim())
1732  {
1733  case 1:
1734  n_dofs = feo->fe<1>()->n_dofs();
1735  break;
1736  case 2:
1737  n_dofs = feo->fe<2>()->n_dofs();
1738  break;
1739  case 3:
1740  n_dofs = feo->fe<3>()->n_dofs();
1741  break;
1742  }
1743 
1744  unsigned int dof_indices[n_dofs];
1745  feo->dh()->get_dof_indices(elem, dof_indices);
1746 
1747  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1748  {
1749  double old_average = 0;
1750  for (unsigned int j=0; j<n_dofs; ++j)
1751  old_average += ls[sbi]->get_solution_array()[dof_indices[j]-feo->dh()->distr()->begin()];
1752  old_average /= n_dofs;
1753 
1754  for (unsigned int j=0; j<n_dofs; ++j)
1755  ls[sbi]->get_solution_array()[dof_indices[j]-feo->dh()->distr()->begin()] += solution_elem_[sbi][i_cell] - old_average;
1756  }
1757  }
1758  }
1759  // update mass_vec for the case that mass matrix changes in next time step
1760  for (unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1761  MatMult(*(ls_dt[sbi]->get_matrix()), ls[sbi]->get_solution(), mass_vec[sbi]);
1762 }
1763 
1764 template<class Model>
1766 {
1767  return Model::mesh_->get_row_4_el();
1768 }
1769 
1770 
1771 
1772 
1773 
1774 
1775 
1776 
1778 template class TransportDG<HeatTransferModel>;
1779 
1780 
1781 
1782 
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
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:50
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:99
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:95
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:189
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
Mat * mass_matrix
The mass matrix.
#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:496
const Vec & get_solution(unsigned int sbi)
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
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:348
EquationOutput output_fields
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
vector< double > sorption_sources
Temporary values of source increments.
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:87
#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
Mat * stiffness_matrix
The stiffness matrix.
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:484
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
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)
Vec * rhs
Vector of right hand side.
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
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
const double epsilon
Definition: mathfce.h:23
FieldCommon & description(const string &description)
Definition: field_common.hh:93
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
virtual MultiFieldValue::return_type value(const Point &p, const ElementAccessor< spacedim > &elm) const
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:59
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:225
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:157
Vec * mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time)...
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.
FieldCommon & name(const string &name)
Definition: field_common.hh:86
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:31
unsigned int el_idx() const
Definition: side_impl.hh:81
void calculate_instant_balance()
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
Definition: linsys.hh:221
Record type proxy class.
Definition: type_record.hh:171
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
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
~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)
Sets up some parameters of the DG method for two sides of an edge.
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
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename MultiFieldValue::return_type > &value_list) const
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