Flow123d  jenkins-Flow123d-windows32-release-multijob-163
transport_dg.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
27  * @author Jan Stebel
28  */
29 
30 #include "system/sys_profiler.hh"
32 
33 #include "io/output_time.hh"
35 #include "fem/mapping_p1.hh"
36 #include "fem/fe_values.hh"
37 #include "fem/fe_p.hh"
38 #include "fem/fe_rt.hh"
39 #include "fields/field_fe.hh"
40 #include "flow/darcy_flow_mh.hh"
41 #include "la/linsys_PETSC.hh"
44 #include "transport/heat_model.hh"
45 
46 #include "fields/generic_field.hh"
47 
48 using namespace Input::Type;
49 
50 template<class Model>
52  = Selection("DG_variant", "Type of penalty term.")
53  .add_value(non_symmetric, "non-symmetric", "non-symmetric weighted interior penalty DG method")
54  .add_value(incomplete, "incomplete", "incomplete weighted interior penalty DG method")
55  .add_value(symmetric, "symmetric", "symmetric weighted interior penalty DG method");
56 
57 template<class Model>
59  Selection("TransportDG_BC_Type", "Types of boundary condition supported by the transport DG model (solute transport or heat transfer).")
60  .add_value(none, "none", "Homogeneous Neumann boundary condition. Zero flux")
61  .add_value(dirichlet, "dirichlet",
62  "Dirichlet boundary condition."
63  //"Specify the pressure head through the 'bc_pressure' field "
64  //"or the piezometric head through the 'bc_piezo_head' field."
65  )
66  .add_value(neumann, "neumann", "Neumann boundary condition. Prescribe water outflow by the 'bc_flux' field.")
67  .add_value(robin, "robin", "Robin boundary condition. Water outflow equal to $\\sigma (h - h^R)$. "
68  //"Specify the transition coefficient by 'bc_sigma' and the reference pressure head or pieaozmetric head "
69  //"through 'bc_pressure' and 'bc_piezo_head' respectively."
70  )
71  .add_value(inflow, "inflow", "Prescribes the concentration in the inflow water on the inflow part of the boundary.");
72 
73 template<class Model>
75  Model::ModelEqData::get_output_selection_input_type("DG", "DG solver")
76  .copy_values(EqData().make_output_field_selection("").close())
77  .close();
78 
79 template<class Model>
81  = Model::get_input_type("DG", "DG solver")
82  .declare_key("solver", LinSys_PETSC::input_type, Default::obligatory(),
83  "Linear solver for MH problem.")
84  .declare_key("input_fields", Array(TransportDG<Model>::EqData().make_field_descriptor_type(std::string(Model::ModelEqData::name()) + "_DG")), IT::Default::obligatory(), "")
85  .declare_key("dg_variant", TransportDG<Model>::dg_variant_selection_input_type, Default("non-symmetric"),
86  "Variant of interior penalty discontinuous Galerkin method.")
87  .declare_key("dg_order", Integer(0,3), Default("1"),
88  "Polynomial order for finite element in DG method (order 0 is suitable if there is no diffusion/dispersion).")
89  .declare_key("output_fields", Array(EqData::output_selection),
90  Default(Model::ModelEqData::default_output_field()),
91  "List of fields to write to output file.");
92 
93 
94 
95 
96 FEObjects::FEObjects(Mesh *mesh_, unsigned int fe_order)
97 {
98  unsigned int q_order;
99 
100  switch (fe_order)
101  {
102  case 0:
103  q_order = 0;
104  fe1_ = new FE_P_disc<0,1,3>;
105  fe2_ = new FE_P_disc<0,2,3>;
106  fe3_ = new FE_P_disc<0,3,3>;
107  break;
108 
109  case 1:
110  q_order = 2;
111  fe1_ = new FE_P_disc<1,1,3>;
112  fe2_ = new FE_P_disc<1,2,3>;
113  fe3_ = new FE_P_disc<1,3,3>;
114  break;
115 
116  case 2:
117  q_order = 4;
118  fe1_ = new FE_P_disc<2,1,3>;
119  fe2_ = new FE_P_disc<2,2,3>;
120  fe3_ = new FE_P_disc<2,3,3>;
121  break;
122 
123  case 3:
124  q_order = 6;
125  fe1_ = new FE_P_disc<3,1,3>;
126  fe2_ = new FE_P_disc<3,2,3>;
127  fe3_ = new FE_P_disc<3,3,3>;
128  break;
129 
130  default:
131  q_order=0;
132  xprintf(PrgErr, "Unsupported polynomial order %d for finite elements in TransportDG ", fe_order);
133  break;
134  }
135 
136  fe_rt1_ = new FE_RT0<1,3>;
137  fe_rt2_ = new FE_RT0<2,3>;
138  fe_rt3_ = new FE_RT0<3,3>;
139 
140  q0_ = new QGauss<0>(q_order);
141  q1_ = new QGauss<1>(q_order);
142  q2_ = new QGauss<2>(q_order);
143  q3_ = new QGauss<3>(q_order);
144 
145  map0_ = new MappingP1<0,3>;
146  map1_ = new MappingP1<1,3>;
147  map2_ = new MappingP1<2,3>;
148  map3_ = new MappingP1<3,3>;
149 
150  dh_ = new DOFHandlerMultiDim(*mesh_);
151 
152  dh_->distribute_dofs(*fe1_, *fe2_, *fe3_);
153 }
154 
155 
157 {
158  delete fe1_;
159  delete fe2_;
160  delete fe3_;
161  delete fe_rt1_;
162  delete fe_rt2_;
163  delete fe_rt3_;
164  delete q0_;
165  delete q1_;
166  delete q2_;
167  delete q3_;
168  delete map0_;
169  delete map1_;
170  delete map2_;
171  delete map3_;
172  delete dh_;
173 }
174 
175 template<> FiniteElement<0,3> *FEObjects::fe<0>() { return 0; }
176 template<> FiniteElement<1,3> *FEObjects::fe<1>() { return fe1_; }
177 template<> FiniteElement<2,3> *FEObjects::fe<2>() { return fe2_; }
178 template<> FiniteElement<3,3> *FEObjects::fe<3>() { return fe3_; }
179 
180 template<> FiniteElement<0,3> *FEObjects::fe_rt<0>() { return 0; }
181 template<> FiniteElement<1,3> *FEObjects::fe_rt<1>() { return fe_rt1_; }
182 template<> FiniteElement<2,3> *FEObjects::fe_rt<2>() { return fe_rt2_; }
183 template<> FiniteElement<3,3> *FEObjects::fe_rt<3>() { return fe_rt3_; }
184 
185 template<> Quadrature<0> *FEObjects::q<0>() { return q0_; }
186 template<> Quadrature<1> *FEObjects::q<1>() { return q1_; }
187 template<> Quadrature<2> *FEObjects::q<2>() { return q2_; }
188 template<> Quadrature<3> *FEObjects::q<3>() { return q3_; }
189 
190 template<> Mapping<0,3> *FEObjects::mapping<0>() { return map0_; }
191 template<> Mapping<1,3> *FEObjects::mapping<1>() { return map1_; }
192 template<> Mapping<2,3> *FEObjects::mapping<2>() { return map2_; }
193 template<> Mapping<3,3> *FEObjects::mapping<3>() { return map3_; }
194 
196 
197 
198 template<class Model>
199 TransportDG<Model>::EqData::EqData() : Model::ModelEqData()
200 {
201  *this+=fracture_sigma
202  .name("fracture_sigma")
203  .description(
204  "Coefficient of diffusive transfer through fractures (for each substance).")
206  .input_default("1.0")
208 
209  *this+=dg_penalty
210  .name("dg_penalty")
211  .description(
212  "Penalty parameter influencing the discontinuity of the solution (for each substance). "
213  "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
215  .input_default("1.0")
217 
218  *this+=bc_type
219  .name("bc_type")
220  .description(
221  "Boundary condition type, possible values: inflow, dirichlet, neumann, robin.")
223  .input_default("\"inflow\"")
225  .flags_add(FieldFlag::in_rhs & FieldFlag::in_main_matrix);
226 
227  *this+=bc_flux
228  .name("bc_flux")
229  .description("Flux in Neumann boundary condition.")
230  .units( UnitSI().kg().m().s(-1).md() )
231  .input_default("0.0")
233  *this+=bc_robin_sigma
234  .name("bc_robin_sigma")
235  .description("Conductivity coefficient in Robin boundary condition.")
236  .units( UnitSI().m(4).s(-1).md() )
237  .input_default("0.0")
238  .flags_add(FieldFlag::in_rhs & FieldFlag::in_main_matrix);
239 
240  *this += region_id.name("region_id")
243 
244  // add all input fields to the output list
245 
246 }
247 
248 template<class Model>
250  : TransportBase(init_mesh, in_rec),
251  mass_matrix(0),
252  allocation_done(false)
253 {
254  // Can not use name() + "constructor" here, since START_TIMER only accepts const char *
255  // due to constexpr optimization.
256  START_TIMER(Model::ModelEqData::name());
257  // Check that Model is derived from AdvectionDiffusionModel.
258  static_assert(std::is_base_of<AdvectionDiffusionModel, Model>::value, "");
259 
260  this->eq_data_ = &data_;
261 
262  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"));
263 
264 
265  // Read names of transported substances.
266  // TODO: Substances should be held in TransportOperatorSplitting only.
267  // Class TransportDG requires only names of components,
268  // and it may have no sense for Model to define Substances
269  // (e.g. if Model represents heat transfer). This should be
270  // resolved when transport classes are refactored so that DG method
271  // can be combined with reactions under operator splitting.
272  Model::set_components(substances_, in_rec);
274 
275  // Set up physical parameters.
276  data_.set_mesh(init_mesh);
277  data_.set_components(substances_.names());
278  data_.set_input_list( in_rec.val<Input::Array>("input_fields") );
279  data_.set_limit_side(LimitSide::right);
281 
282 
283  // DG variant and order
284  dg_variant = in_rec.val<DGVariant>("dg_variant");
285  dg_order = in_rec.val<unsigned int>("dg_order");
286 
287  // DG stabilization parameters on boundary edges
288  gamma.resize(n_subst_);
289  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
290  gamma[sbi].resize(mesh_->boundary_.size());
291 
292  // create finite element structures and distribute DOFs
293  feo = new FEObjects(mesh_, dg_order);
294  DBGMSG("TDG: solution size %d\n", feo->dh()->n_global_dofs());
295 
296  // Resize coefficient arrays
297  int qsize = max(feo->q<0>()->size(), max(feo->q<1>()->size(), max(feo->q<2>()->size(), feo->q<3>()->size())));
298  int max_edg_sides = max(mesh_->max_edge_sides(1), max(mesh_->max_edge_sides(2), mesh_->max_edge_sides(3)));
299  mm_coef.resize(qsize);
300  ad_coef.resize(n_subst_);
301  dif_coef.resize(n_subst_);
302  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
303  {
304  ad_coef[sbi].resize(qsize);
305  dif_coef[sbi].resize(qsize);
306  }
307  ad_coef_edg.resize(max_edg_sides);
308  dif_coef_edg.resize(max_edg_sides);
309  for (int sd=0; sd<max_edg_sides; sd++)
310  {
311  ad_coef_edg[sd].resize(n_subst_);
312  dif_coef_edg[sd].resize(n_subst_);
313  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
314  {
315  ad_coef_edg[sd][sbi].resize(qsize);
316  dif_coef_edg[sd][sbi].resize(qsize);
317  }
318  }
319 
320  // register output fields
321  output_rec = in_rec.val<Input::Record>("output_stream");
322  output_vec.resize(n_subst_);
323  output_solution.resize(n_subst_);
324  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
325  {
326  // for each substance we allocate output array and vector
327  output_solution[sbi] = new double[feo->dh()->n_global_dofs()];
328  VecCreateSeqWithArray(PETSC_COMM_SELF, 1, feo->dh()->n_global_dofs(), output_solution[sbi], &output_vec[sbi]);
329  }
330  data_.output_field.set_components(substances_.names());
331  data_.output_field.set_mesh(*mesh_);
332  data_.output_type(OutputTime::CORNER_DATA);
333 
334  data_.output_field.set_up_components();
335  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
336  {
337  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
338  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar> > output_field_ptr(new FieldFE<3, FieldValue<3>::Scalar>);
339  output_field_ptr->set_fe_data(feo->dh(), feo->mapping<1>(), feo->mapping<2>(), feo->mapping<3>(), &output_vec[sbi]);
340  data_.output_field[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
341  }
342  data_.set_limit_side(LimitSide::left);
344  output_stream->add_admissible_field_names(in_rec.val<Input::Array>("output_fields"));
345 
346  // set time marks for writing the output
348 
349  // allocate matrix and vector structures
350  ls = new LinSys*[n_subst_];
351  ls_dt = new LinSys_PETSC(feo->dh()->distr());
352  ( (LinSys_PETSC *)ls_dt )->set_from_input( in_rec.val<Input::Record>("solver") );
353  for (unsigned int sbi = 0; sbi < n_subst_; sbi++) {
354  ls[sbi] = new LinSys_PETSC(feo->dh()->distr());
355  ( (LinSys_PETSC *)ls[sbi] )->set_from_input( in_rec.val<Input::Record>("solver") );
356  ls[sbi]->set_solution(NULL);
357  }
358  stiffness_matrix = new Mat[n_subst_];
359  rhs = new Vec[n_subst_];
360  mass_vec = new Vec[n_subst_];
361 
362 
363  // initialization of balance object
365  if (it->val<bool>("balance_on"))
366  {
367  balance_ = boost::make_shared<Balance>(Model::balance_prefix(), mesh_, feo->dh()->el_ds(), feo->dh()->get_el_4_loc(), *it);
368 
369  // a not very nice workaround for model with single solution component with no name
370  if (typeid(Model) == typeid(HeatTransferModel))
371  subst_idx = {balance_->add_quantity("energy")};
372  else
373  subst_idx = balance_->add_quantities(substances_.names());
374 
375  balance_->allocate(feo->dh()->distr()->lsize(),
376  max(feo->fe<1>()->n_dofs(), max(feo->fe<2>()->n_dofs(), feo->fe<3>()->n_dofs())));
377  }
378 
379 }
380 
381 template<class Model>
383 {
384  delete time_;
385  delete ls_dt;
386 
387  if (feo->dh()->el_ds()->myp() == 0)
388  {
389  for (unsigned int i=0; i<n_subst_; i++)
390  {
391  VecDestroy(&output_vec[i]);
392  delete[] output_solution[i];
393  }
394  }
395 
396  for (unsigned int i=0; i<n_subst_; i++)
397  {
398  delete ls[i];
399  MatDestroy(&stiffness_matrix[i]);
400  VecDestroy(&rhs[i]);
401  VecDestroy(&mass_vec[i]);
402  }
403  delete[] ls;
404  delete[] stiffness_matrix;
405  delete[] rhs;
406  delete[] mass_vec;
407  delete feo;
408 
409  gamma.clear();
410  delete output_stream;
411 }
412 
413 
414 template<class Model>
416 {
417  IS is;
418  VecScatter output_scatter;
419  int idx[] = { 0 };
420  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
421  {
422  // gather solution to output_vec[sbi]
423  ISCreateBlock(PETSC_COMM_SELF, ls[sbi]->size(), 1, idx, PETSC_COPY_VALUES, &is);
424  VecScatterCreate(ls[sbi]->get_solution(), is, output_vec[sbi], PETSC_NULL, &output_scatter);
425  VecScatterBegin(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
426  VecScatterEnd(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
427  VecScatterDestroy(&(output_scatter));
428  ISDestroy(&(is));
429  }
430 }
431 
432 
433 
434 template<class Model>
436 {
437  START_TIMER(Model::ModelEqData::name());
438  data_.mark_input_times(time_->equation_fixed_mark_type());
439  data_.set_time(time_->step());
440 
441 
442  // set initial conditions
443  set_initial_condition();
444  for (unsigned int sbi = 0; sbi < n_subst_; sbi++)
445  ( (LinSys_PETSC *)ls[sbi] )->set_initial_guess_nonzero();
446 
447  // during preallocation we assemble the matrices and vectors required for mass balance
448  if (balance_ != nullptr)
449  {
450  balance_->units(Model::balance_units());
451 
452  if (!allocation_done) preallocate();
453 
454  for (unsigned int sbi=0; sbi<n_subst_; ++sbi)
455  {
456  balance_->calculate_mass(subst_idx[sbi], ls[sbi]->get_solution());
457  balance_->calculate_source(subst_idx[sbi], ls[sbi]->get_solution());
458  balance_->calculate_flux(subst_idx[sbi], ls[sbi]->get_solution());
459  }
460  }
461 
462  output_data();
463 }
464 
465 
466 template<class Model>
468 {
469  // preallocate mass matrix
470  ls_dt->start_allocation();
471  assemble_mass_matrix();
472  mass_matrix = NULL;
473 
474  // preallocate system matrix
475  for (unsigned int i=0; i<n_subst_; i++)
476  {
477  ls[i]->start_allocation();
478  stiffness_matrix[i] = NULL;
479  rhs[i] = NULL;
480  }
481  assemble_stiffness_matrix();
482  set_sources();
483  set_boundary_conditions();
484 
485  allocation_done = true;
486 }
487 
488 
489 
490 template<class Model>
492 {
493  START_TIMER("DG-ONE STEP");
494 
495  time_->next_time();
496  time_->view("TDG");
497 
498  START_TIMER("data reinit");
499  data_.set_time(time_->step());
500  END_TIMER("data reinit");
501 
502  // check first time assembly - needs preallocation
503  if (!allocation_done) preallocate();
504 
505  // assemble mass matrix
506  if (mass_matrix == NULL || data_.subset(FieldFlag::in_time_term).changed() )
507  {
508  ls_dt->start_add_assembly();
509  ls_dt->mat_zero_entries();
510  assemble_mass_matrix();
511  ls_dt->finish_assembly();
512 
513  // construct mass_vec for initial time
514  if (mass_matrix == NULL)
515  {
516  for (unsigned int i=0; i<n_subst_; i++)
517  {
518  VecDuplicate(ls[i]->get_solution(), &mass_vec[i]);
519  MatMult(*(ls_dt->get_matrix()), ls[i]->get_solution(), mass_vec[i]);
520  }
521  }
522 
523  mass_matrix = *(ls_dt->get_matrix());
524  }
525 
526  // assemble stiffness matrix
527  if (stiffness_matrix[0] == NULL
528  || data_.subset(FieldFlag::in_main_matrix).changed()
529  || Model::flux_changed)
530  {
531  // new fluxes can change the location of Neumann boundary,
532  // thus stiffness matrix must be reassembled
533  for (unsigned int i=0; i<n_subst_; i++)
534  {
535  ls[i]->start_add_assembly();
536  ls[i]->mat_zero_entries();
537  }
538  assemble_stiffness_matrix();
539  for (unsigned int i=0; i<n_subst_; i++)
540  {
541  ls[i]->finish_assembly();
542 
543  if (stiffness_matrix[i] == NULL)
544  MatConvert(*( ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
545  else
546  MatCopy(*( ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
547  }
548  }
549 
550  // assemble right hand side (due to sources and boundary conditions)
551  if (rhs[0] == NULL
552  || data_.subset(FieldFlag::in_rhs).changed()
553  || Model::flux_changed)
554  {
555  for (unsigned int i=0; i<n_subst_; i++)
556  {
557  ls[i]->start_add_assembly();
558  ls[i]->rhs_zero_entries();
559  }
560  set_sources();
561  set_boundary_conditions();
562  for (unsigned int i=0; i<n_subst_; i++)
563  {
564  ls[i]->finish_assembly();
565 
566  VecDuplicate(*( ls[i]->get_rhs() ), &rhs[i]);
567  VecCopy(*( ls[i]->get_rhs() ), rhs[i]);
568  }
569  }
570 
571  Model::flux_changed = false;
572 
573 
574  /* Apply backward Euler time integration.
575  *
576  * Denoting A the stiffness matrix and M the mass matrix, the algebraic system at the k-th time level reads
577  *
578  * (1/dt M + A)u^k = f + 1/dt M.u^{k-1}
579  *
580  * Hence we modify at each time level the right hand side:
581  *
582  * f^k = f + 1/dt M u^{k-1},
583  *
584  * where f stands for the term stemming from the force and boundary conditions.
585  * Accordingly, we set
586  *
587  * A^k = A + 1/dt M.
588  *
589  */
590  Mat m;
591  START_TIMER("solve");
592  for (unsigned int i=0; i<n_subst_; i++)
593  {
594  MatConvert(stiffness_matrix[i], MATSAME, MAT_INITIAL_MATRIX, &m);
595  MatAXPY(m, 1./time_->dt(), mass_matrix, SUBSET_NONZERO_PATTERN);
596  ls[i]->set_matrix(m, DIFFERENT_NONZERO_PATTERN);
597  Vec w;
598  VecDuplicate(rhs[i], &w);
599  VecWAXPY(w, 1./time_->dt(), mass_vec[i], rhs[i]);
600  ls[i]->set_rhs(w);
601 
602  VecDestroy(&w);
603  MatDestroy(&m);
604 
605  ls[i]->solve();
606 
607  MatMult(*(ls_dt->get_matrix()), ls[i]->get_solution(), mass_vec[i]);
608  }
609  END_TIMER("solve");
610 
611  if (balance_ != nullptr)
612  {
613  for (unsigned int sbi=0; sbi<n_subst_; ++sbi)
614  {
615  balance_->calculate_mass(subst_idx[sbi], ls[sbi]->get_solution());
616  balance_->calculate_source(subst_idx[sbi], ls[sbi]->get_solution());
617  balance_->calculate_flux(subst_idx[sbi], ls[sbi]->get_solution());
618  if (balance_->cumulative())
619  {
620  balance_->calculate_cumulative_sources(subst_idx[sbi], ls[sbi]->get_solution(), time_->dt());
621  balance_->calculate_cumulative_fluxes(subst_idx[sbi], ls[sbi]->get_solution(), time_->dt());
622  }
623  }
624  }
625 
626  END_TIMER("DG-ONE STEP");
627 }
628 
629 
630 template<class Model>
632 {
633  // So far the velocity_vector contains zeros, so we ignore it.
634  // Instead we use the value Side.flux.
635 
636  mh_dh = &dh;
637  Model::flux_changed = true;
638 
639 }
640 
641 
642 template<class Model>
644 {
645  if (!time_->is_current( time_->marks().type_output() )) return;
646 
647  START_TIMER("DG-OUTPUT");
648 
649  // gather the solution from all processors
650  output_vector_gather();
651  data_.subset(FieldFlag::allow_output).set_time( time_->step());
652  data_.output(output_stream);
653  output_stream->write_time_frame();
654 
655  if (balance_ != nullptr)
656  balance_->output(time_->t());
657 
658  END_TIMER("DG-OUTPUT");
659 }
660 
661 
662 template<class Model>
664 {
665  START_TIMER("assemble_mass");
666  if (balance_ != nullptr)
667  balance_->start_mass_assembly(subst_idx);
668  assemble_mass_matrix<1>();
669  assemble_mass_matrix<2>();
670  assemble_mass_matrix<3>();
671  if (balance_ != nullptr)
672  balance_->finish_mass_assembly(subst_idx);
673  END_TIMER("assemble_mass");
674 }
675 
676 
677 template<class Model> template<unsigned int dim>
679 {
680  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(), update_values | update_JxW_values | update_quadrature_points);
681  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
682  vector<int> dof_indices(ndofs);
683  PetscScalar local_mass_matrix[ndofs*ndofs];
684  vector<PetscScalar> local_mass_balance_vector(ndofs);
685 
686  // assemble integral over elements
687  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
688  {
689  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
690  if (cell->dim() != dim) continue;
691 
692  fe_values.reinit(cell);
693  feo->dh()->get_dof_indices(cell, (unsigned int*)&(dof_indices[0]));
694  ElementAccessor<3> ele_acc = cell->element_accessor();
695 
696  Model::compute_mass_matrix_coefficient(fe_values.point_list(), ele_acc, mm_coef);
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]*fe_values.shape_value(j,k)*fe_values.shape_value(i,k)*fe_values.JxW(k);
706  }
707 
708  if (balance_ != nullptr)
709  {
710  local_mass_balance_vector[i] = 0;
711  for (unsigned int k=0; k<qsize; k++)
712  local_mass_balance_vector[i] += mm_coef[k]*fe_values.shape_value(i,k)*fe_values.JxW(k);
713  }
714  }
715 
716  if (balance_ != nullptr)
717  for (unsigned int sbi=0; sbi<n_subst_; ++sbi)
718  balance_->add_mass_matrix_values(subst_idx[sbi], ele_acc.region().bulk_idx(), dof_indices, local_mass_balance_vector);
719 
720  ls_dt->mat_set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), local_mass_matrix);
721  }
722 }
723 
724 
725 
726 
727 template<class Model>
729 {
730  START_TIMER("assemble_stiffness");
731  START_TIMER("assemble_volume_integrals");
732  assemble_volume_integrals<1>();
733  assemble_volume_integrals<2>();
734  assemble_volume_integrals<3>();
735  END_TIMER("assemble_volume_integrals");
736 
737  START_TIMER("assemble_fluxes_boundary");
738  assemble_fluxes_boundary<1>();
739  assemble_fluxes_boundary<2>();
740  assemble_fluxes_boundary<3>();
741  END_TIMER("assemble_fluxes_boundary");
742 
743  START_TIMER("assemble_fluxes_elem_elem");
744  assemble_fluxes_element_element<1>();
745  assemble_fluxes_element_element<2>();
746  assemble_fluxes_element_element<3>();
747  END_TIMER("assemble_fluxes_elem_elem");
748 
749  START_TIMER("assemble_fluxes_elem_side");
750  assemble_fluxes_element_side<1>();
751  assemble_fluxes_element_side<2>();
752  assemble_fluxes_element_side<3>();
753  END_TIMER("assemble_fluxes_elem_side");
754  END_TIMER("assemble_stiffness");
755 }
756 
757 
758 
759 template<class Model>
760 template<unsigned int dim>
762 {
763  FEValues<dim,3> fv_rt(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe_rt<dim>(),
765  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
767  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
768  unsigned int dof_indices[ndofs];
769  vector<arma::vec3> velocity(qsize);
770  vector<arma::vec> sources_sigma(qsize, arma::vec(n_substances()));
771  PetscScalar local_matrix[ndofs*ndofs];
772 
773  // assemble integral over elements
774  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
775  {
776  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
777  if (cell->dim() != dim) continue;
778 
779  fe_values.reinit(cell);
780  fv_rt.reinit(cell);
781  ElementAccessor<3> ele_acc = cell->element_accessor();
782  feo->dh()->get_dof_indices(cell, dof_indices);
783 
784  calculate_velocity(cell, velocity, fv_rt);
785  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), velocity, ele_acc, ad_coef, dif_coef);
786  Model::compute_sources_sigma(fe_values.point_list(), ele_acc, sources_sigma);
787 
788  // assemble the local stiffness matrix
789  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
790  {
791  for (unsigned int i=0; i<ndofs; i++)
792  for (unsigned int j=0; j<ndofs; j++)
793  local_matrix[i*ndofs+j] = 0;
794 
795  for (unsigned int k=0; k<qsize; k++)
796  {
797  for (unsigned int i=0; i<ndofs; i++)
798  {
799  arma::vec3 Kt_grad_i = dif_coef[sbi][k].t()*fe_values.shape_grad(i,k);
800  double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.shape_grad(i,k));
801 
802  for (unsigned int j=0; j<ndofs; j++)
803  local_matrix[i*ndofs+j] += (arma::dot(Kt_grad_i, fe_values.shape_grad(j,k))
804  -fe_values.shape_value(j,k)*ad_dot_grad_i
805  +sources_sigma[k][sbi]*fe_values.shape_value(j,k)*fe_values.shape_value(i,k))*fe_values.JxW(k);
806  }
807  }
808  ls[sbi]->mat_set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, local_matrix);
809  }
810  }
811 }
812 
813 
814 template<class Model>
816 {
817  START_TIMER("assemble_sources");
818  if (balance_ != nullptr)
819  balance_->start_source_assembly(subst_idx);
820  set_sources<1>();
821  set_sources<2>();
822  set_sources<3>();
823  if (balance_ != nullptr)
824  balance_->finish_source_assembly(subst_idx);
825  END_TIMER("assemble_sources");
826 }
827 
828 template<class Model>
829 template<unsigned int dim>
831 {
832  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
834  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
835  vector<arma::vec> sources_conc(qsize, arma::vec(n_substances())),
836  sources_density(qsize, arma::vec(n_substances())),
837  sources_sigma(qsize, arma::vec(n_substances()));
838  vector<int> dof_indices(ndofs);
839  PetscScalar local_rhs[ndofs];
840  vector<PetscScalar> local_source_balance_vector(ndofs), local_source_balance_rhs(ndofs);
841  double source;
842 
843  // assemble integral over elements
844  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
845  {
846  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
847  if (cell->dim() != dim) continue;
848 
849  fe_values.reinit(cell);
850  feo->dh()->get_dof_indices(cell, (unsigned int *)&(dof_indices[0]));
851 
852  Model::compute_source_coefficients(fe_values.point_list(), cell->element_accessor(), sources_conc, sources_density, sources_sigma);
853 
854  // assemble the local stiffness matrix
855  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
856  {
857  fill_n(local_rhs, ndofs, 0);
858  local_source_balance_vector.assign(ndofs, 0);
859  local_source_balance_rhs.assign(ndofs, 0);
860 
861  // compute sources
862  for (unsigned int k=0; k<qsize; k++)
863  {
864  source = (sources_density[k][sbi] + sources_conc[k][sbi]*sources_sigma[k][sbi])*fe_values.JxW(k);
865 
866  for (unsigned int i=0; i<ndofs; i++)
867  {
868  local_rhs[i] += source*fe_values.shape_value(i,k);
869 
870  if (balance_ != nullptr)
871  {
872  local_source_balance_vector[i] -= sources_sigma[k][sbi]*fe_values.shape_value(i,k)*fe_values.JxW(k);
873  local_source_balance_rhs[i] += source*fe_values.shape_value(i,k);
874  }
875  }
876  }
877  ls[sbi]->rhs_set_values(ndofs, &(dof_indices[0]), local_rhs);
878 
879  if (balance_ != nullptr)
880  {
881  balance_->add_source_matrix_values(subst_idx[sbi], cell->region().bulk_idx(), dof_indices, local_source_balance_vector);
882  balance_->add_source_rhs_values(subst_idx[sbi], cell->region().bulk_idx(), dof_indices, local_source_balance_rhs);
883  }
884  }
885  }
886 }
887 
888 
889 
890 template<class Model>
891 template<unsigned int dim>
893 {
894  vector<FESideValues<dim,3>*> fe_values;
895  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
896  update_values);
897  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size(),
898  n_max_sides = ad_coef_edg.size();
899  vector<unsigned int*> side_dof_indices;
900  PetscScalar local_matrix[ndofs*ndofs];
901  vector<vector<arma::vec3> > side_velocity(n_max_sides);
902  vector<arma::vec> dg_penalty(n_max_sides);
903  double gamma_l, omega[2], transport_flux;
904 
905  for (unsigned int sid=0; sid<n_max_sides; sid++)
906  {
907  side_dof_indices.push_back(new unsigned int[ndofs]);
908  fe_values.push_back(new FESideValues<dim,3>(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
910  }
911 
912  // assemble integral over sides
913  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
914  {
915  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
916  if (edg->n_sides < 2 || edg->side(0)->element()->dim() != dim) continue;
917 
918  for (int sid=0; sid<edg->n_sides; sid++)
919  {
920  typename DOFHandlerBase::CellIterator cell = edg->side(sid)->element();
921  ElementAccessor<3> ele_acc = cell->element_accessor();
922  feo->dh()->get_dof_indices(cell, side_dof_indices[sid]);
923  fe_values[sid]->reinit(cell, edg->side(sid)->el_idx());
924  fsv_rt.reinit(cell, edg->side(sid)->el_idx());
925  calculate_velocity(cell, side_velocity[sid], fsv_rt);
926  Model::compute_advection_diffusion_coefficients(fe_values[sid]->point_list(), side_velocity[sid], ele_acc, ad_coef_edg[sid], dif_coef_edg[sid]);
927  dg_penalty[sid] = data_.dg_penalty.value(cell->centre(), ele_acc);
928  }
929 
930  // fluxes and penalty
931  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
932  {
933  vector<double> fluxes(edg->n_sides);
934  for (int sid=0; sid<edg->n_sides; sid++)
935  {
936  fluxes[sid] = 0;
937  for (unsigned int k=0; k<qsize; k++)
938  fluxes[sid] += arma::dot(ad_coef_edg[sid][sbi][k], fe_values[sid]->normal_vector(k))*fe_values[sid]->JxW(k);
939  fluxes[sid] /= edg->side(sid)->measure();
940  }
941 
942  for (int s1=0; s1<edg->n_sides; s1++)
943  {
944  for (int s2=s1+1; s2<edg->n_sides; s2++)
945  {
946  ASSERT(edg->side(s1)->valid(), "Invalid side of edge.");
947  if (!feo->dh()->el_is_local(edg->side(s1)->element().index())
948  && !feo->dh()->el_is_local(edg->side(s2)->element().index())) continue;
949 
950  arma::vec3 nv = fe_values[s1]->normal_vector(0);
951 
952  // set up the parameters for DG method
953  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);
954 
955  int sd[2];
956  sd[0] = s1;
957  sd[1] = s2;
958 
959 #define AVERAGE(i,k,side_id) (fe_values[sd[side_id]]->shape_value(i,k)*0.5)
960 #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])
961 #define JUMP(i,k,side_id) ((side_id==0?1:-1)*fe_values[sd[side_id]]->shape_value(i,k))
962 
963  // For selected pair of elements:
964  for (int n=0; n<2; n++)
965  {
966  if (!feo->dh()->el_is_local(edg->side(sd[n])->element().index())) continue;
967 
968  for (int m=0; m<2; m++)
969  {
970  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
971  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
972  local_matrix[i*fe_values[sd[m]]->n_dofs()+j] = 0;
973 
974  for (unsigned int k=0; k<qsize; k++)
975  {
976  double flux_times_JxW = transport_flux*fe_values[0]->JxW(k);
977  double gamma_times_JxW = gamma_l*fe_values[0]->JxW(k);
978 
979  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
980  {
981  double flux_JxW_jump_i = flux_times_JxW*JUMP(i,k,n);
982  double gamma_JxW_jump_i = gamma_times_JxW*JUMP(i,k,n);
983  double JxW_jump_i = fe_values[0]->JxW(k)*JUMP(i,k,n);
984  double JxW_var_wavg_i = fe_values[0]->JxW(k)*WAVERAGE(i,k,n)*dg_variant;
985 
986  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
987  {
988  int index = i*fe_values[sd[m]]->n_dofs()+j;
989 
990  // flux due to transport (applied on interior edges) (average times jump)
991  local_matrix[index] += flux_JxW_jump_i*AVERAGE(j,k,m);
992 
993  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
994  local_matrix[index] += gamma_JxW_jump_i*JUMP(j,k,m);
995 
996  // terms due to diffusion
997  local_matrix[index] -= WAVERAGE(j,k,m)*JxW_jump_i;
998  local_matrix[index] -= JUMP(j,k,m)*JxW_var_wavg_i;
999  }
1000  }
1001  }
1002  ls[sbi]->mat_set_values(fe_values[sd[n]]->n_dofs(), (int *)side_dof_indices[sd[n]], fe_values[sd[m]]->n_dofs(), (int *)side_dof_indices[sd[m]], local_matrix);
1003  }
1004  }
1005 #undef AVERAGE
1006 #undef WAVERAGE
1007 #undef JUMP
1008  }
1009  }
1010  }
1011  }
1012 
1013  for (unsigned int i=0; i<n_max_sides; i++)
1014  {
1015  delete fe_values[i];
1016  delete[] side_dof_indices[i];
1017  }
1018 }
1019 
1020 
1021 template<class Model>
1022 template<unsigned int dim>
1024 {
1025  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1027  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1028  update_values);
1029  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1030  unsigned int side_dof_indices[ndofs];
1031  PetscScalar local_matrix[ndofs*ndofs];
1032  vector<arma::vec3> side_velocity;
1033  vector<arma::vec> robin_sigma(qsize, arma::vec(n_substances()));
1034  arma::vec dg_penalty;
1035  double gamma_l;
1036 
1037  // assemble boundary integral
1038  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1039  {
1040  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
1041  if (edg->n_sides > 1) continue;
1042  // check spatial dimension
1043  if (edg->side(0)->dim() != dim-1) continue;
1044  // skip edges lying not on the boundary
1045  if (edg->side(0)->cond() == NULL) continue;
1046 
1047  SideIter side = edg->side(0);
1048  typename DOFHandlerBase::CellIterator cell = side->element();
1049  ElementAccessor<3> ele_acc = cell->element_accessor();
1050  feo->dh()->get_dof_indices(cell, side_dof_indices);
1051  fe_values_side.reinit(cell, side->el_idx());
1052  fsv_rt.reinit(cell, side->el_idx());
1053 
1054  calculate_velocity(cell, side_velocity, fsv_rt);
1055  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), side_velocity, ele_acc, ad_coef, dif_coef);
1056  dg_penalty = data_.dg_penalty.value(cell->centre(), ele_acc);
1057  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), side->cond()->element_accessor());
1058  data_.bc_robin_sigma.value_list(fe_values_side.point_list(), side->cond()->element_accessor(), robin_sigma);
1059 
1060  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1061  {
1062  for (unsigned int i=0; i<ndofs; i++)
1063  for (unsigned int j=0; j<ndofs; j++)
1064  local_matrix[i*ndofs+j] = 0;
1065 
1066  // On Neumann boundaries we have only term from integrating by parts the advective term,
1067  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
1068  double side_flux = 0;
1069  for (unsigned int k=0; k<qsize; k++)
1070  side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
1071  double transport_flux = side_flux/side->measure();
1072 
1073  if (bc_type[sbi] == EqData::dirichlet)
1074  {
1075  // set up the parameters for DG method
1076  set_DG_parameters_boundary(side, qsize, dif_coef[sbi], transport_flux, fe_values_side.normal_vector(0), dg_penalty[sbi], gamma_l);
1077  gamma[sbi][side->cond_idx()] = gamma_l;
1078  transport_flux += gamma_l;
1079  }
1080 
1081  // fluxes and penalty
1082  for (unsigned int k=0; k<qsize; k++)
1083  {
1084  double flux_times_JxW;
1085  if (bc_type[sbi] == EqData::robin)
1086  flux_times_JxW = (transport_flux + robin_sigma[k][sbi])*fe_values_side.JxW(k);
1087  else if (bc_type[sbi] == EqData::inflow && side_flux < 0)
1088  flux_times_JxW = 0;
1089  else
1090  flux_times_JxW = transport_flux*fe_values_side.JxW(k);
1091 
1092  for (unsigned int i=0; i<ndofs; i++)
1093  {
1094  for (unsigned int j=0; j<ndofs; j++)
1095  {
1096  // flux due to advection and penalty
1097  local_matrix[i*ndofs+j] += flux_times_JxW*fe_values_side.shape_value(i,k)*fe_values_side.shape_value(j,k);
1098 
1099  // flux due to diffusion (only on dirichlet and inflow boundary)
1100  if (bc_type[sbi] == EqData::dirichlet)
1101  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)
1102  + 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
1103  )*fe_values_side.JxW(k);
1104  }
1105  }
1106  }
1107 
1108  ls[sbi]->mat_set_values(ndofs, (int *)side_dof_indices, ndofs, (int *)side_dof_indices, local_matrix);
1109  }
1110  }
1111 }
1112 
1113 
1114 template<class Model>
1115 template<unsigned int dim>
1117 {
1118 
1119  if (dim == 1) return;
1120  FEValues<dim-1,3> fe_values_vb(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe<dim-1>(),
1122  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1124  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1125  update_values);
1126  FEValues<dim-1,3> fv_rt(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe_rt<dim-1>(),
1127  update_values);
1128 
1129  vector<FEValuesSpaceBase<3>*> fv_sb(2);
1130  const unsigned int ndofs = feo->fe<dim>()->n_dofs(); // number of local dofs
1131  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
1132  unsigned int side_dof_indices[2*ndofs], n_dofs[2];
1133  vector<arma::vec3> velocity_higher, velocity_lower;
1134  vector<arma::vec> frac_sigma(qsize, arma::vec(n_substances()));
1135  vector<double> csection_lower(qsize), csection_higher(qsize), mm_coef_lower(qsize), mm_coef_higher(qsize);
1136  PetscScalar local_matrix[4*ndofs*ndofs];
1137  double comm_flux[2][2];
1138 
1139  // index 0 = element with lower dimension,
1140  // index 1 = side of element with higher dimension
1141  fv_sb[0] = &fe_values_vb;
1142  fv_sb[1] = &fe_values_side;
1143 
1144  // assemble integral over sides
1145  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
1146  {
1147  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
1148  // skip neighbours of different dimension
1149  if (nb->element()->dim() != dim-1) continue;
1150 
1151  typename DOFHandlerBase::CellIterator cell_sub = mesh_->element.full_iter(nb->element());
1152  feo->dh()->get_dof_indices(cell_sub, side_dof_indices);
1153  fe_values_vb.reinit(cell_sub);
1154  n_dofs[0] = fv_sb[0]->n_dofs();
1155 
1156  typename DOFHandlerBase::CellIterator cell = nb->side()->element();
1157  feo->dh()->get_dof_indices(cell, side_dof_indices+n_dofs[0]);
1158  fe_values_side.reinit(cell, nb->side()->el_idx());
1159  n_dofs[1] = fv_sb[1]->n_dofs();
1160 
1161  // Element id's for testing if they belong to local partition.
1162  int element_id[2];
1163  element_id[0] = cell_sub.index();
1164  element_id[1] = cell.index();
1165 
1166  fsv_rt.reinit(cell, nb->side()->el_idx());
1167  fv_rt.reinit(cell_sub);
1168  calculate_velocity(cell, velocity_higher, fsv_rt);
1169  calculate_velocity(cell_sub, velocity_lower, fv_rt);
1170  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_lower, cell_sub->element_accessor(), ad_coef_edg[0], dif_coef_edg[0]);
1171  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_higher, cell->element_accessor(), ad_coef_edg[1], dif_coef_edg[1]);
1172  Model::compute_mass_matrix_coefficient(fe_values_vb.point_list(), cell_sub->element_accessor(), mm_coef_lower);
1173  Model::compute_mass_matrix_coefficient(fe_values_vb.point_list(), cell->element_accessor(), mm_coef_higher);
1174  data_.cross_section.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), csection_lower);
1175  data_.cross_section.value_list(fe_values_vb.point_list(), cell->element_accessor(), csection_higher);
1176  data_.fracture_sigma.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), frac_sigma);
1177 
1178  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1179  {
1180  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
1181  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
1182  local_matrix[i*(n_dofs[0]+n_dofs[1])+j] = 0;
1183 
1184  // set transmission conditions
1185  for (unsigned int k=0; k<qsize; k++)
1186  {
1187  /* The communication flux has two parts:
1188  * - "diffusive" term containing sigma
1189  * - "advective" term representing usual upwind
1190  *
1191  * The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
1192  * than b and A in the manual.
1193  * In calculation of sigma there appears one more csection_lower in the denominator.
1194  */
1195  double sigma = frac_sigma[k][sbi]*arma::dot(dif_coef_edg[0][sbi][k]*fe_values_side.normal_vector(k),fe_values_side.normal_vector(k))*
1196  2*csection_higher[k]*csection_higher[k]/(csection_lower[k]*csection_lower[k]);
1197 
1198  // Since mm_coef_* contains cross section, we have to divide by it.
1199  double transport_flux = arma::dot(ad_coef_edg[1][sbi][k], fe_values_side.normal_vector(k));
1200  double por_lower_over_higher = mm_coef_lower[k]*csection_higher[k]/(mm_coef_higher[k]*csection_lower[k]);
1201 
1202  comm_flux[0][0] = (sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1203  comm_flux[0][1] = -(sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1204  comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1205  comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1206 
1207  for (int n=0; n<2; n++)
1208  {
1209  if (!feo->dh()->el_is_local(element_id[n])) continue;
1210 
1211  for (unsigned int i=0; i<n_dofs[n]; i++)
1212  for (int m=0; m<2; m++)
1213  for (unsigned int j=0; j<n_dofs[m]; j++)
1214  local_matrix[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
1215  comm_flux[m][n]*fv_sb[m]->shape_value(j,k)*fv_sb[n]->shape_value(i,k);
1216  }
1217  }
1218  ls[sbi]->mat_set_values(n_dofs[0]+n_dofs[1], (int *)side_dof_indices, n_dofs[0]+n_dofs[1], (int *)side_dof_indices, local_matrix);
1219  }
1220  }
1221 
1222 }
1223 
1224 
1225 
1226 
1227 
1228 
1229 template<class Model>
1231 {
1232  START_TIMER("assemble_bc");
1233  if (balance_ != nullptr)
1234  balance_->start_flux_assembly(subst_idx);
1235  set_boundary_conditions<1>();
1236  set_boundary_conditions<2>();
1237  set_boundary_conditions<3>();
1238  if (balance_ != nullptr)
1239  balance_->finish_flux_assembly(subst_idx);
1240  END_TIMER("assemble_bc");
1241 }
1242 
1243 
1244 template<class Model>
1245 template<unsigned int dim>
1247 {
1248  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1250  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1251  update_values);
1252  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1253  vector<int> side_dof_indices(ndofs);
1254  unsigned int loc_b=0;
1255  double local_rhs[ndofs];
1256  vector<PetscScalar> local_flux_balance_vector(ndofs);
1257  PetscScalar local_flux_balance_rhs;
1258  vector<arma::vec> bc_values(qsize, arma::vec(n_substances())),
1259  bc_fluxes(qsize, arma::vec(n_substances())),
1260  bc_sigma(qsize, arma::vec(n_substances()));
1261  vector<arma::vec3> velocity;
1262 
1263  for (unsigned int loc_el = 0; loc_el < feo->dh()->el_ds()->lsize(); loc_el++)
1264  {
1265  ElementFullIter elm = mesh_->element(feo->dh()->el_index(loc_el));
1266  if (elm->boundary_idx_ == nullptr) continue;
1267 
1268  FOR_ELEMENT_SIDES(elm,si)
1269  {
1270  Edge *edg = elm->side(si)->edge();
1271  if (edg->n_sides > 1) continue;
1272  // skip edges lying not on the boundary
1273  if (edg->side(0)->cond() == NULL) continue;
1274 
1275  if (edg->side(0)->dim() != dim-1)
1276  {
1277  if (edg->side(0)->cond() != nullptr) ++loc_b;
1278  continue;
1279  }
1280 
1281  SideIter side = edg->side(0);
1282  typename DOFHandlerBase::CellIterator cell = mesh().element.full_iter(side->element());
1283  ElementAccessor<3> ele_acc = side->cond()->element_accessor();
1284 
1285  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), ele_acc);
1286 
1287  fe_values_side.reinit(cell, side->el_idx());
1288  fsv_rt.reinit(cell, side->el_idx());
1289  calculate_velocity(cell, velocity, fsv_rt);
1290 
1291  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), velocity, side->element()->element_accessor(), ad_coef, dif_coef);
1292  Model::compute_dirichlet_bc(fe_values_side.point_list(), ele_acc, bc_values);
1293  data_.bc_flux.value_list(fe_values_side.point_list(), ele_acc, bc_fluxes);
1294  data_.bc_robin_sigma.value_list(fe_values_side.point_list(), ele_acc, bc_sigma);
1295 
1296  feo->dh()->get_dof_indices(cell, (unsigned int *)&(side_dof_indices[0]));
1297 
1298  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1299  {
1300  fill_n(local_rhs, ndofs, 0);
1301  local_flux_balance_vector.assign(ndofs, 0);
1302  local_flux_balance_rhs = 0;
1303 
1304  double side_flux = 0;
1305  for (unsigned int k=0; k<qsize; k++)
1306  side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
1307  double transport_flux = side_flux/side->measure();
1308 
1309  for (unsigned int k=0; k<qsize; k++)
1310  {
1311  double bc_term = 0;
1312  arma::vec3 bc_grad;
1313  bc_grad.zeros();
1314  if (bc_type[sbi] == EqData::inflow && side_flux < 0)
1315  {
1316  bc_term = -transport_flux*bc_values[k][sbi]*fe_values_side.JxW(k);
1317  }
1318  else if (bc_type[sbi] == EqData::dirichlet)
1319  {
1320  bc_term = gamma[sbi][side->cond_idx()]*bc_values[k][sbi]*fe_values_side.JxW(k);
1321  bc_grad = -bc_values[k][sbi]*fe_values_side.JxW(k)*dg_variant*(arma::trans(dif_coef[sbi][k])*fe_values_side.normal_vector(k));
1322  }
1323  else if (bc_type[sbi] == EqData::neumann)
1324  {
1325  bc_term = -bc_fluxes[k][sbi]*fe_values_side.JxW(k);
1326  }
1327  else if (bc_type[sbi] == EqData::robin)
1328  {
1329  bc_term = bc_sigma[k][sbi]*bc_values[k][sbi]*fe_values_side.JxW(k);
1330  }
1331 
1332  for (unsigned int i=0; i<ndofs; i++)
1333  {
1334  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k)
1335  + arma::dot(bc_grad,fe_values_side.shape_grad(i,k));
1336 
1337  if (balance_ != nullptr)
1338  {
1339  if (bc_type[sbi] == EqData::dirichlet)
1340  {
1341  local_flux_balance_vector[i] += (arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.shape_value(i,k)
1342  - arma::dot(dif_coef[sbi][k]*fe_values_side.shape_grad(i,k),fe_values_side.normal_vector(k))
1343  + gamma[sbi][side->cond_idx()]*fe_values_side.shape_value(i,k))*fe_values_side.JxW(k);
1344  local_flux_balance_rhs -= (time_->tlevel() == 0?0:1)*bc_term*fe_values_side.shape_value(i,k);
1345  }
1346  else if (bc_type[sbi] == EqData::inflow && side_flux < 0)
1347  {
1348  local_flux_balance_rhs -= bc_term*fe_values_side.shape_value(i,k);
1349  }
1350  else if (bc_type[sbi] == EqData::neumann)
1351  {
1352  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);
1353  local_flux_balance_rhs -= bc_term*fe_values_side.shape_value(i,k);
1354  }
1355  else if (bc_type[sbi] == EqData::robin)
1356  {
1357  local_flux_balance_vector[i] += (arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k)) + bc_sigma[k][sbi])*fe_values_side.JxW(k)*fe_values_side.shape_value(i,k);
1358  local_flux_balance_rhs -= bc_term*fe_values_side.shape_value(i,k);
1359  }
1360  else
1361  {
1362  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);
1363  }
1364  }
1365  }
1366  }
1367  ls[sbi]->rhs_set_values(ndofs, &(side_dof_indices[0]), local_rhs);
1368 
1369  if (balance_ != nullptr)
1370  {
1371  balance_->add_flux_matrix_values(subst_idx[sbi], loc_b, side_dof_indices, local_flux_balance_vector);
1372  balance_->add_flux_vec_value(subst_idx[sbi], loc_b, local_flux_balance_rhs);
1373  }
1374  }
1375  ++loc_b;
1376  }
1377  }
1378 }
1379 
1380 
1381 
1382 // TODO: The detection of side number from SideIter
1383 // in TransportDG::calculate_velocity()
1384 // should be done with help of class RefElement. This however requires
1385 // that the MH_DofHandler uses the node/side ordering defined in
1386 // the respective RefElement.
1387 template<class Model>
1388 template<unsigned int dim>
1390 {
1391  std::map<const Node*, int> node_nums;
1392  for (unsigned int i=0; i<cell->n_nodes(); i++)
1393  node_nums[cell->node[i]] = i;
1394 
1395  velocity.resize(fv.n_points());
1396 
1397  for (unsigned int k=0; k<fv.n_points(); k++)
1398  {
1399  velocity[k].zeros();
1400  for (unsigned int sid=0; sid<cell->n_sides(); sid++)
1401  {
1402  if (cell->side(sid)->dim() != dim-1) continue;
1403  int num = dim*(dim+1)/2;
1404  for (unsigned int i=0; i<cell->side(sid)->n_nodes(); i++)
1405  num -= node_nums[cell->side(sid)->node(i)];
1406  velocity[k] += fv.shape_vector(num,k) * mh_dh->side_flux( *(cell->side(sid)) );
1407  }
1408  }
1409 }
1410 
1411 
1412 
1413 
1414 
1415 template<class Model>
1417  const int s1,
1418  const int s2,
1419  const int K_size,
1420  const vector<arma::mat33> &K1,
1421  const vector<arma::mat33> &K2,
1422  const vector<double> &fluxes,
1423  const arma::vec3 &normal_vector,
1424  const double alpha1,
1425  const double alpha2,
1426  double &gamma,
1427  double *omega,
1428  double &transport_flux)
1429 {
1430  double delta[2];
1431  double h = 0;
1432  double local_alpha = 0;
1433 
1434  ASSERT(edg.side(s1)->valid(), "Invalid side of an edge.");
1435  SideIter s = edg.side(s1);
1436 
1437  // calculate the side diameter
1438  if (s->dim() == 0)
1439  {
1440  h = 1;
1441  }
1442  else
1443  {
1444  for (unsigned int i=0; i<s->n_nodes(); i++)
1445  for (unsigned int j=i+1; j<s->n_nodes(); j++)
1446  h = max(h, s->node(i)->distance(*s->node(j)));
1447  }
1448 
1449  // calculate the total in- and out-flux through the edge
1450  double pflux = 0, nflux = 0;
1451  for (int i=0; i<edg.n_sides; i++)
1452  {
1453  if (fluxes[i] > 0)
1454  pflux += fluxes[i];
1455  else
1456  nflux += fluxes[i];
1457  }
1458 
1459  // calculate the flux from s1 to s2
1460  if (fluxes[s2] > 0 && fluxes[s1] < 0 && s1 < s2)
1461  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
1462  else if (fluxes[s2] < 0 && fluxes[s1] > 0 && s1 < s2)
1463  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
1464  else if (s1==s2)
1465  transport_flux = fluxes[s1];
1466  else
1467  transport_flux = 0;
1468 
1469  gamma = 0.5*fabs(transport_flux);
1470 
1471 
1472  // determine local DG penalty
1473  local_alpha = max(alpha1, alpha2);
1474 
1475  if (s1 == s2)
1476  {
1477  omega[0] = 1;
1478 
1479  // delta is set to the average value of Kn.n on the side
1480  delta[0] = 0;
1481  for (int k=0; k<K_size; k++)
1482  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1483  delta[0] /= K_size;
1484 
1485  gamma += local_alpha/h*(delta[0]);
1486  }
1487  else
1488  {
1489  delta[0] = 0;
1490  delta[1] = 0;
1491  for (int k=0; k<K_size; k++)
1492  {
1493  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1494  delta[1] += dot(K2[k]*normal_vector,normal_vector);
1495  }
1496  delta[0] /= K_size;
1497  delta[1] /= K_size;
1498 
1499  double delta_sum = delta[0] + delta[1];
1500 
1501  if (delta_sum > numeric_limits<double>::epsilon())
1502  {
1503  omega[0] = delta[1]/delta_sum;
1504  omega[1] = delta[0]/delta_sum;
1505  gamma += local_alpha/h*(delta[0]*delta[1]/delta_sum);
1506  }
1507  else
1508  for (int i=0; i<2; i++) omega[i] = 0;
1509  }
1510 }
1511 
1512 
1513 
1514 
1515 
1516 
1517 template<class Model>
1519  const int K_size,
1520  const vector<arma::mat33> &K,
1521  const double flux,
1522  const arma::vec3 &normal_vector,
1523  const double alpha,
1524  double &gamma)
1525 {
1526  double delta = 0, h = 0;
1527 
1528  // calculate the side diameter
1529  if (side->dim() == 0)
1530  {
1531  h = 1;
1532  }
1533  else
1534  {
1535  for (unsigned int i=0; i<side->n_nodes(); i++)
1536  for (unsigned int j=i+1; j<side->n_nodes(); j++)
1537  h = max(h, side->node(i)->distance( *side->node(j) ));
1538  }
1539 
1540  // delta is set to the average value of Kn.n on the side
1541  for (int k=0; k<K_size; k++)
1542  delta += dot(K[k]*normal_vector,normal_vector);
1543  delta /= K_size;
1544 
1545  gamma = 0.5*fabs(flux) + alpha/h*delta;
1546 }
1547 
1548 
1549 
1550 
1551 
1552 template<class Model>
1554 {
1555  START_TIMER("set_init_cond");
1556  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1557  ls[sbi]->start_allocation();
1558  prepare_initial_condition<1>();
1559  prepare_initial_condition<2>();
1560  prepare_initial_condition<3>();
1561 
1562  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1563  ls[sbi]->start_add_assembly();
1564  prepare_initial_condition<1>();
1565  prepare_initial_condition<2>();
1566  prepare_initial_condition<3>();
1567 
1568  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1569  {
1570  ls[sbi]->finish_assembly();
1571  ls[sbi]->solve();
1572  }
1573  END_TIMER("set_init_cond");
1574 }
1575 
1576 template<class Model>
1577 template<unsigned int dim>
1579 {
1580  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1582  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1583  unsigned int dof_indices[ndofs];
1584  double matrix[ndofs*ndofs], rhs[ndofs];
1585  std::vector<arma::vec> init_values(qsize);
1586 
1587  for (unsigned int k=0; k<qsize; k++)
1588  init_values[k].resize(n_subst_);
1589 
1590  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
1591  {
1592  typename DOFHandlerBase::CellIterator elem = mesh_->element(feo->dh()->el_index(i_cell));
1593  if (elem->dim() != dim) continue;
1594 
1595  ElementAccessor<3> ele_acc = elem->element_accessor();
1596  feo->dh()->get_dof_indices(elem, dof_indices);
1597  fe_values.reinit(elem);
1598 
1599  Model::compute_init_cond(fe_values.point_list(), ele_acc, init_values);
1600 
1601  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1602  {
1603  for (unsigned int i=0; i<ndofs; i++)
1604  {
1605  rhs[i] = 0;
1606  for (unsigned int j=0; j<ndofs; j++)
1607  matrix[i*ndofs+j] = 0;
1608  }
1609 
1610  for (unsigned int k=0; k<qsize; k++)
1611  {
1612  double rhs_term = init_values[k](sbi)*fe_values.JxW(k);
1613 
1614  for (unsigned int i=0; i<ndofs; i++)
1615  {
1616  for (unsigned int j=0; j<ndofs; j++)
1617  matrix[i*ndofs+j] += fe_values.shape_value(i,k)*fe_values.shape_value(j,k)*fe_values.JxW(k);
1618 
1619  rhs[i] += fe_values.shape_value(i,k)*rhs_term;
1620  }
1621  }
1622  ls[sbi]->set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, matrix, rhs);
1623  }
1624  }
1625 }
1626 
1627 
1628 template<class Model>
1629 void TransportDG<Model>::calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance)
1630 {
1631  calc_fluxes<1>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1632  calc_fluxes<2>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1633  calc_fluxes<3>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1634 }
1635 
1636 template<class Model>
1637 template<unsigned int dim>
1638 void TransportDG<Model>::calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance)
1639 {
1640  FESideValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1642  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1643  update_values);
1644  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1645  unsigned int dof_indices[ndofs];
1646  vector<arma::vec3> side_velocity(qsize);
1647  vector<arma::vec> bc_values(qsize, arma::vec(n_subst_)), bc_fluxes(qsize, arma::vec(n_subst_)), bc_sigma(qsize, arma::vec(n_subst_));
1648  arma::vec3 conc_grad;
1649 
1650  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1651  {
1652  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
1653  if (edg->n_sides > 1) continue;
1654  if (edg->side(0)->dim() != dim-1) continue;
1655  // skip edges lying not on the boundary
1656  if (edg->side(0)->cond() == NULL) continue;
1657 
1658  SideIter side = edg->side(0);
1660  ElementAccessor<3> ele_acc = cell->element_accessor();
1661  Region r = side->cond()->element_accessor().region();
1662  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d, edg: % d\n", side->cond()->bc_ele_idx_, side->cond()->edge_idx_);
1663  unsigned int bc_region_idx = r.boundary_idx();
1664 
1665  fe_values.reinit(cell, side->el_idx());
1666  fsv_rt.reinit(cell, side->el_idx());
1667  feo->dh()->get_dof_indices(cell, dof_indices);
1668  calculate_velocity(cell, side_velocity, fsv_rt);
1669  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), side_velocity, ele_acc, ad_coef, dif_coef);
1670  Model::compute_dirichlet_bc(fe_values.point_list(), side->cond()->element_accessor(), bc_values);
1671  data_.bc_flux.value_list(fe_values.point_list(), side->cond()->element_accessor(), bc_fluxes);
1672  data_.bc_robin_sigma.value_list(fe_values.point_list(), side->cond()->element_accessor(), bc_sigma);
1673  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), side->cond()->element_accessor());
1674 
1675  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1676  {
1677  double mass_flux = 0;
1678  double water_flux = 0;
1679  for (unsigned int k=0; k<qsize; k++)
1680  water_flux += arma::dot(ad_coef[sbi][k], fe_values.normal_vector(k))*fe_values.JxW(k);
1681  water_flux /= side->measure();
1682 
1683  for (unsigned int k=0; k<qsize; k++)
1684  {
1685  double conc = 0;
1686  conc_grad.zeros();
1687  for (unsigned int i=0; i<ndofs; i++)
1688  {
1689  conc += fe_values.shape_value(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1690  conc_grad += fe_values.shape_grad(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1691  }
1692 
1693  // flux due to advection
1694  if (bc_type[sbi] == EqData::inflow && water_flux < 0)
1695  mass_flux += water_flux*bc_values[k][sbi]*fe_values.JxW(k);
1696  else
1697  mass_flux += water_flux*conc*fe_values.JxW(k);
1698 
1699  if (bc_type[sbi] == EqData::dirichlet)
1700  {
1701  // flux due to diffusion
1702  mass_flux -= arma::dot(dif_coef[sbi][k]*conc_grad,fe_values.normal_vector(k))*fe_values.JxW(k);
1703 
1704  // the penalty term has to be added otherwise the mass balance will not hold
1705  mass_flux -= gamma[sbi][side->cond_idx()]*(bc_values[k][sbi] - conc)*fe_values.JxW(k);
1706  }
1707  else if (bc_type[sbi] == EqData::neumann)
1708  {
1709  // flux due to Neumann b.c.
1710  mass_flux += bc_fluxes[k][sbi]*fe_values.JxW(k);
1711  }
1712  else if (bc_type[sbi] == EqData::robin)
1713  {
1714  // flux due to Robin b.c.
1715  mass_flux += bc_sigma[k][sbi]*(conc - bc_values[k][sbi])*fe_values.JxW(k);
1716  }
1717 
1718  }
1719  bcd_balance[sbi][bc_region_idx] += mass_flux;
1720  if (mass_flux >= 0) bcd_plus_balance[sbi][bc_region_idx] += mass_flux;
1721  else bcd_minus_balance[sbi][bc_region_idx] += mass_flux;
1722  }
1723  }
1724 
1725 }
1726 
1727 template<class Model>
1729 {
1730  calc_elem_sources<1>(mass, src_balance);
1731  calc_elem_sources<2>(mass, src_balance);
1732  calc_elem_sources<3>(mass, src_balance);
1733 }
1734 
1735 template<class Model>
1736 template<unsigned int dim>
1738 {
1739  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1741  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1742  unsigned int dof_indices[ndofs];
1743  vector<arma::vec> sources_conc(qsize, arma::vec(n_substances())),
1744  sources_density(qsize, arma::vec(n_substances())),
1745  sources_sigma(qsize, arma::vec(n_substances()));
1746  double mass_sum, sources_sum, conc, conc_diff;
1747 
1748  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
1749  {
1750  typename DOFHandlerBase::CellIterator elem = mesh_->element(feo->dh()->el_index(i_cell));
1751  if (elem->dim() != dim) continue;
1752 
1753  Region r = elem->element_accessor().region();
1754  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d\n", elem.index());
1755  unsigned int region_idx = r.bulk_idx();
1756 
1757  ElementAccessor<3> ele_acc = elem->element_accessor();
1758  fe_values.reinit(elem);
1759  feo->dh()->get_dof_indices(elem, dof_indices);
1760 
1761  Model::compute_mass_matrix_coefficient(fe_values.point_list(), ele_acc, mm_coef);
1762  Model::compute_source_coefficients(fe_values.point_list(), ele_acc, sources_conc, sources_density, sources_sigma);
1763 
1764  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1765  {
1766  mass_sum = 0;
1767  sources_sum = 0;
1768 
1769  for (unsigned int k=0; k<qsize; k++)
1770  {
1771  conc = 0;
1772  for (unsigned int i=0; i<ndofs; i++)
1773  conc += fe_values.shape_value(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1774 
1775  mass_sum += mm_coef[k]*conc*fe_values.JxW(k);
1776 
1777  conc_diff = sources_conc[k][sbi] - conc;
1778  sources_sum += (sources_density[k][sbi] + conc_diff*sources_sigma[k][sbi])*fe_values.JxW(k);
1779  }
1780 
1781  mass[sbi][region_idx] += mass_sum;
1782  src_balance[sbi][region_idx] += sources_sum;
1783  }
1784  }
1785 }
1786 
1787 
1788 
1789 
1790 
1792 template class TransportDG<HeatTransferModel>;
1793 
1794 
1795 
1796 
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
Shape function values.
Definition: update_flags.hh:98
FieldSet * eq_data_
Definition: equation.hh:227
Distribution * el_ds() const
Returns the distribution of number of element to local processes.
Definition: dofhandler.hh:306
#define DBGMSG(...)
Definition: global_defs.h:196
void assemble_fluxes_element_element()
Assembles the fluxes between elements of the same dimension.
void set_sources()
Assembles the right hand side due to volume sources.
void set_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
Transformed quadrature weight for cell sides.
FiniteElement< dim, 3 > * fe()
vector< double > mm_coef
Mass matrix coefficients.
Accessor to input data conforming to declared Array.
Definition: accessors.hh:558
void assemble_fluxes_element_side()
Assembles the fluxes between elements of different dimensions.
Definition: system.hh:72
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:34
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
Input::Record output_rec
Record with output specification.
static constexpr Mask allow_output
The field can output. Is part of generated output selection. (default on)
Definition: field_flag.hh:27
int tlevel() const
const std::vector< std::string > & names()
Definition: substance.hh:97
Transport with dispersion implemented using discontinuous Galerkin method.
Field< 3, FieldValue< 3 >::Integer > region_id
void add_admissible_field_names(const Input::Array &in_array)
Registers names of output fields that can be written using this stream.
Definition: output_time.cc:165
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:39
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
double distance(const Node &n2) const
Definition: nodes.hh:98
Base class for FEValues and FESideValues.
Definition: fe_values.hh:43
void assemble_mass_matrix()
Assembles the mass matrix.
Boundary * cond() const
Definition: side_impl.hh:59
OutputTime * output_stream
int * get_el_4_loc() const
Definition: dofhandler.hh:308
void next_time()
Proceed to the next time according to current estimated time step.
void update_solution() override
Computes the solution in one time instant.
TimeMark::Type type_output()
Definition: time_marks.hh:203
static Default obligatory()
Definition: type_record.hh:89
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...
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:361
const MH_DofHandler * mh_dh
Definition: mesh.h:109
Fields computed from the mesh data.
Iterator< Ret > find(const string &key) const
void set_initial_condition()
Sets the initial condition.
Class FEValues calculates finite element data on the actual cells such as shape function values...
#define FOR_ELEMENT_SIDES(i, j)
Definition: elements.h:151
int index() const
Definition: sys_vector.hh:88
bool is_current(const TimeMark::Type &mask) const
int n_sides
Definition: edges.h:48
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
static Input::Type::Record input_type
Definition: linsys_PETSC.hh:46
Definition: edges.h:38
bool valid() const
Definition: side_impl.hh:75
const RegionDB & region_db() const
Definition: mesh.h:155
void set_from_input(const Input::Record in_rec)
#define AVERAGE(i, k, side_id)
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
const TimeStep & step(int index=-1) const
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:355
double t() const
vector< Boundary > boundary_
Definition: mesh.h:209
Basic time management functionality for unsteady (and steady) solvers (class Equation).
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
unsigned int size() const
Definition: substance.hh:99
static TimeMarks & marks()
const double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Definition: fe_values.cc:178
Class for declaration of inputs sequences.
Definition: type_base.hh:239
Specification of transport model interface.
void view(const char *name="") const
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
Symmetric Gauss-Legendre quadrature formulae on simplices.
#define WAVERAGE(i, k, side_id)
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
Definition: region.hh:75
unsigned int dim() const
TransportDG(Mesh &init_mesh, const Input::Record &in_rec)
Constructor.
vector< double * > output_solution
Array for storing the output solution data.
vector< vector< vector< arma::vec3 > > > ad_coef_edg
Advection coefficients on edges.
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:213
Transformed quadrature points.
void output_data()
Postprocesses the solution and writes to output file.
#define ASSERT(...)
Definition: global_defs.h:121
static OutputTime * create_output_stream(const Input::Record &in_rec)
This method delete all object instances of class OutputTime stored in output_streams vector...
Definition: output_time.cc:141
void preallocate()
unsigned int dim() const
Definition: side_impl.hh:26
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
FEObjects * feo
Finite element objects.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:40
Mesh & mesh()
Definition: equation.hh:171
static constexpr Mask in_time_term
A field is part of time term of the equation.
Definition: field_flag.hh:32
FieldCommon & input_default(const string &input_default)
Definition: field_common.hh:92
TimeMark::Type equation_fixed_mark_type() const
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:44
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:42
SideIter side()
Selection & add_value(const int value, const std::string &key, const std::string &description="")
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
FieldCommon & input_selection(const Input::Type::Selection *element_selection)
const Ret val(const string &key) const
static auto region_id(Mesh &mesh) -> IndexField
#define xprintf(...)
Definition: system.hh:100
Mat mass_matrix
The mass matrix.
#define START_TIMER(tag)
Starts a timer with specified tag.
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:251
Mat * stiffness_matrix
The stiffness matrix.
Mesh * mesh_
Definition: equation.hh:218
Element * element()
Definition: boundaries.cc:47
Field< 3, FieldValue< 3 >::Vector > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
boost::shared_ptr< Balance > balance_
(new) object for calculation and writing the mass balance to file.
unsigned int dg_order
Polynomial order of finite elements.
void output_vector_gather()
const vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.cc:190
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:36
unsigned int n_subst_
Number of transported substances.
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:64
Shape function gradients.
void mark_output_times(const TimeGovernor &tg)
Definition: output_time.cc:182
Vec * rhs
Vector of right hand side.
void set_solution(double *sol_array)
Definition: linsys.hh:264
unsigned int edge_idx_
Definition: boundaries.h:87
BCField< 3, FieldValue< 3 >::Vector > bc_flux
Flux in Neumann or Robin b.c.
static Input::Type::Selection bc_type_selection
ElementIter element()
double measure() const
Definition: sides.cc:42
const 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:154
DOFHandlerMultiDim * dh()
Region region() const
Definition: accessors.hh:88
Normal vectors.
unsigned int bc_ele_idx_
Definition: boundaries.h:88
LinSys * ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
Discontinuous Galerkin method for equation of transport with dispersion.
const double epsilon
Definition: mathfce.h:6
void output_data() override
Write computed fields.
vector< Neighbour > vb_neighbours_
Definition: mesh.h:235
double side_flux(const Side &side) const
temporary replacement for DofHandler accessor, flux through given side
FieldCommon & description(const string &description)
Definition: field_common.hh:80
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:64
std::vector< std::vector< double > > gamma
Penalty parameters.
unsigned int max_edge_sides(unsigned int dim) const
Definition: mesh.h:182
Discontinuous Galerkin method for equation of transport with dispersion.
SubstanceList substances_
Transported substances.
ElementFullIter element() const
Definition: side_impl.hh:41
void calc_fluxes(vector< vector< double > > &bcd_balance, vector< vector< double > > &bcd_plus_balance, vector< vector< double > > &bcd_minus_balance)
Calculates flux through boundary of each region.
double dt() const
arma::vec3 centre() const
Definition: elements.cc:132
Definition: system.hh:72
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:312
BCField< 3, FieldValue< 3 >::EnumVector > bc_type
Type of boundary condition (see also BC_Type)
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:244
Vec * mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time)...
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:216
const arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.cc:196
FiniteElement< dim, spacedim > * fe
The used finite element.
Definition: fe_values.hh:337
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:80
Abstract linear system class.
Definitions of particular quadrature rules on simplices.
FieldCommon & name(const string &name)
Definition: field_common.hh:73
void calc_elem_sources(vector< vector< double > > &mass, vector< vector< double > > &src_balance)
Calculates volume sources for each region.
vector< Vec > output_vec
Vector of solution data.
Calculates finite element data on the actual cell.
Definition: fe_values.hh:370
#define END_TIMER(tag)
Ends a timer with specified tag.
Discontinuous Galerkin method for equation of transport with dispersion.
Quadrature< dim > * q()
virtual void set_velocity_field(const MH_DofHandler &dh)
Updates the velocity field which determines some coefficients of the transport equation.
mixed-hybrid model of linear Darcy flow, possibly unsteady.
unsigned int el_idx() const
Definition: side_impl.hh:70
const 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:148
unsigned int n_substances() override
Returns number of trnasported substances.
Record type proxy class.
Definition: type_record.hh:169
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
bool is_valid() const
Returns false if the region has undefined/invalid value.
Definition: region.hh:67
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Class for representation SI units of Fields.
Definition: unit_si.hh:31
const unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.cc:202
const Node * node(unsigned int i) const
Definition: side_impl.hh:35
const arma::vec::fixed< spacedim > shape_vector(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Definition: fe_values.cc:160
void calculate_velocity(const ElementFullIter &cell, std::vector< arma::vec3 > &velocity, FEValuesBase< dim, 3 > &fv)
Calculates the velocity field on a given dim dimensional cell.
Distribution * distr() const
Definition: dofhandler.hh:90
Mapping< dim, 3 > * mapping()
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:44
~TransportDG()
Destructor.
LinSys ** ls
Linear algebra system for the transport equation.
void set_DG_parameters_edge(const Edge &edg, const int s1, const int s2, const int K_size, const std::vector< arma::mat33 > &K1, const std::vector< arma::mat33 > &K2, const std::vector< double > &fluxes, const arma::vec3 &normal_vector, const double alpha1, const double alpha2, double &gamma, double *omega, double &transport_flux)
Sets up some parameters of the DG method for two sides of an edge.
SideIter side(const unsigned int i) const
Definition: edges.h:43
Field< 3, FieldValue< 3 >::Vector > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
#define JUMP(i, k, side_id)
Template for classes storing finite set of named values.
void set_DG_parameters_boundary(const SideIter side, const int K_size, const std::vector< arma::mat33 > &K, const double flux, const arma::vec3 &normal_vector, const double alpha, double &gamma)
Sets up parameters of the DG method on a given boundary edge.
Definitions of Raviart-Thomas finite elements.
const unsigned int n_global_dofs() const
Getter for the number of all mesh dofs required by the given finite element.
Definition: dofhandler.hh:72
BCField< 3, FieldValue< 3 >::Vector > bc_robin_sigma
Transition coefficient in Robin b.c.
FullIter full_iter(Iter it)
Definition: sys_vector.hh:438
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:55
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)
Definition: transport_dg.cc:96
TimeGovernor * time_
Definition: equation.hh:219
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:205
unsigned int n_nodes() const
Definition: side_impl.hh:22
Transformed quadrature weights.
Calculates finite element data on a side.
Definition: fe_values.hh:415
unsigned int lsize(int proc) const
get local size