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