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