Flow123d  jenkins-Flow123d-windows32-release-multijob-51
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"
33 #include "fem/mapping_p1.hh"
34 #include "fem/fe_values.hh"
35 #include "fem/fe_p.hh"
36 #include "fem/fe_rt.hh"
37 #include "io/output.h"
38 #include "fields/field_fe.hh"
39 #include "flow/darcy_flow_mh.hh"
40 #include "la/linsys_PETSC.hh"
43 #include "transport/heat_model.hh"
44 
45 using namespace Input::Type;
46 
47 template<class Model>
49  = Selection("DG_variant", "Type of penalty term.")
50  .add_value(non_symmetric, "non-symmetric", "non-symmetric weighted interior penalty DG method")
51  .add_value(incomplete, "incomplete", "incomplete weighted interior penalty DG method")
52  .add_value(symmetric, "symmetric", "symmetric weighted interior penalty DG method");
53 
54 template<class Model>
56  Selection("TransportDG_EqData_bc_Type")
57  .add_value(inflow, "inflow", "Dirichlet BC on inflow and homogeneous Neumann BC on outflow.")
58  .add_value(dirichlet, "dirichlet")
59  .add_value(neumann, "neumann")
60  .add_value(robin, "robin");
61 
62 template<class Model>
64  Model::ModelEqData::get_output_selection_input_type("DG", "DG solver")
65  .copy_values(EqData().make_output_field_selection("").close())
66  .close();
67 
68 template<class Model>
70  = Model::get_input_type("DG", "DG solver")
71  .declare_key("solver", LinSys_PETSC::input_type, Default::obligatory(),
72  "Linear solver for MH problem.")
73  .declare_key("input_fields", Array(TransportDG<Model>::EqData().make_field_descriptor_type(Model::ModelEqData::name() + "_DG")), IT::Default::obligatory(), "")
74  .declare_key("dg_variant", TransportDG<Model>::dg_variant_selection_input_type, Default("non-symmetric"),
75  "Variant of interior penalty discontinuous Galerkin method.")
76  .declare_key("dg_order", Integer(0,3), Default("1"),
77  "Polynomial order for finite element in DG method (order 0 is suitable if there is no diffusion/dispersion).")
78  .declare_key("output_fields", Array(EqData::output_selection),
79  Default(Model::ModelEqData::default_output_field()),
80  "List of fields to write to output file.");
81 
82 
83 
84 FEObjects::FEObjects(Mesh *mesh_, unsigned int fe_order)
85 {
86  unsigned int q_order;
87 
88  switch (fe_order)
89  {
90  case 0:
91  q_order = 0;
92  fe1_ = new FE_P_disc<0,1,3>;
93  fe2_ = new FE_P_disc<0,2,3>;
94  fe3_ = new FE_P_disc<0,3,3>;
95  break;
96 
97  case 1:
98  q_order = 2;
99  fe1_ = new FE_P_disc<1,1,3>;
100  fe2_ = new FE_P_disc<1,2,3>;
101  fe3_ = new FE_P_disc<1,3,3>;
102  break;
103 
104  case 2:
105  q_order = 4;
106  fe1_ = new FE_P_disc<2,1,3>;
107  fe2_ = new FE_P_disc<2,2,3>;
108  fe3_ = new FE_P_disc<2,3,3>;
109  break;
110 
111  case 3:
112  q_order = 6;
113  fe1_ = new FE_P_disc<3,1,3>;
114  fe2_ = new FE_P_disc<3,2,3>;
115  fe3_ = new FE_P_disc<3,3,3>;
116  break;
117 
118  default:
119  q_order=0;
120  xprintf(PrgErr, "Unsupported polynomial order %d for finite elements in TransportDG ", fe_order);
121  break;
122  }
123 
124  fe_rt1_ = new FE_RT0<1,3>;
125  fe_rt2_ = new FE_RT0<2,3>;
126  fe_rt3_ = new FE_RT0<3,3>;
127 
128  q0_ = new QGauss<0>(q_order);
129  q1_ = new QGauss<1>(q_order);
130  q2_ = new QGauss<2>(q_order);
131  q3_ = new QGauss<3>(q_order);
132 
133  map0_ = new MappingP1<0,3>;
134  map1_ = new MappingP1<1,3>;
135  map2_ = new MappingP1<2,3>;
136  map3_ = new MappingP1<3,3>;
137 
138  dh_ = new DOFHandlerMultiDim(*mesh_);
139 
140  dh_->distribute_dofs(*fe1_, *fe2_, *fe3_);
141 }
142 
143 
145 {
146  delete fe1_;
147  delete fe2_;
148  delete fe3_;
149  delete fe_rt1_;
150  delete fe_rt2_;
151  delete fe_rt3_;
152  delete q0_;
153  delete q1_;
154  delete q2_;
155  delete q3_;
156  delete map0_;
157  delete map1_;
158  delete map2_;
159  delete map3_;
160  delete dh_;
161 }
162 
163 template<> FiniteElement<0,3> *FEObjects::fe<0>() { return 0; }
164 template<> FiniteElement<1,3> *FEObjects::fe<1>() { return fe1_; }
165 template<> FiniteElement<2,3> *FEObjects::fe<2>() { return fe2_; }
166 template<> FiniteElement<3,3> *FEObjects::fe<3>() { return fe3_; }
167 
168 template<> FiniteElement<0,3> *FEObjects::fe_rt<0>() { return 0; }
169 template<> FiniteElement<1,3> *FEObjects::fe_rt<1>() { return fe_rt1_; }
170 template<> FiniteElement<2,3> *FEObjects::fe_rt<2>() { return fe_rt2_; }
171 template<> FiniteElement<3,3> *FEObjects::fe_rt<3>() { return fe_rt3_; }
172 
173 template<> Quadrature<0> *FEObjects::q<0>() { return q0_; }
174 template<> Quadrature<1> *FEObjects::q<1>() { return q1_; }
175 template<> Quadrature<2> *FEObjects::q<2>() { return q2_; }
176 template<> Quadrature<3> *FEObjects::q<3>() { return q3_; }
177 
178 template<> Mapping<0,3> *FEObjects::mapping<0>() { return map0_; }
179 template<> Mapping<1,3> *FEObjects::mapping<1>() { return map1_; }
180 template<> Mapping<2,3> *FEObjects::mapping<2>() { return map2_; }
181 template<> Mapping<3,3> *FEObjects::mapping<3>() { return map3_; }
182 
184 
185 
186 template<class Model>
187 TransportDG<Model>::EqData::EqData() : Model::ModelEqData()
188 {
189  *this+=fracture_sigma
190  .name("fracture_sigma")
191  .description(
192  "Coefficient of diffusive transfer through fractures (for each substance).")
194  .input_default("1.0")
196 
197  *this+=dg_penalty
198  .name("dg_penalty")
199  .description(
200  "Penalty parameter influencing the discontinuity of the solution (for each substance). "
201  "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
203  .input_default("1.0")
205 
206  *this+=bc_type
207  .name("bc_type")
208  .description(
209  "Boundary condition type, possible values: inflow, dirichlet, neumann, robin.")
211  .input_default("\"inflow\"")
213  .flags_add(FieldFlag::in_rhs & FieldFlag::in_main_matrix);
214 
215  *this+=bc_flux
216  .name("bc_flux")
217  .description("Flux in Neumann boundary condition.")
218  .units( UnitSI().kg().m().s(-1).md() )
219  .input_default("0.0")
221  *this+=bc_robin_sigma
222  .name("bc_robin_sigma")
223  .description("Conductivity coefficient in Robin boundary condition.")
224  .units( UnitSI().m(4).s(-1).md() )
225  .input_default("0.0")
226  .flags_add(FieldFlag::in_rhs & FieldFlag::in_main_matrix);
227 
228  // add all input fields to the output list
229 
230 }
231 
232 template<class Model>
234  : TransportBase(init_mesh, in_rec),
235  mass_matrix(0),
236  allocation_done(false)
237 {
238  // Check that Model is derived from AdvectionDiffusionModel.
239  static_assert(std::is_base_of<AdvectionDiffusionModel, Model>::value, "");
240 
241  this->eq_data_ = &data_;
242 
243  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"));
245 
246  // Read names of transported substances.
247  Model::set_component_names(subst_names_, in_rec);
248  n_subst_ = subst_names_.size();
249 
250  Input::Iterator<Input::Record> it = in_rec.find<Input::Record>("mass_balance");
251  if (it)
252  mass_balance_ = new MassBalance(this, *it);
253 
254  // Set up physical parameters.
255  data_.set_mesh(init_mesh);
256  data_.set_n_components(n_subst_);
257  data_.set_input_list( in_rec.val<Input::Array>("input_fields") );
258  data_.set_limit_side(LimitSide::left);
259 
260 
261  // DG variant and order
262  dg_variant = in_rec.val<DGVariant>("dg_variant");
263  dg_order = in_rec.val<unsigned int>("dg_order");
264 
265  // DG stabilization parameters on boundary edges
266  gamma.resize(n_subst_);
267  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
268  gamma[sbi].resize(mesh_->boundary_.size());
269 
270  // create finite element structures and distribute DOFs
271  feo = new FEObjects(mesh_, dg_order);
272  DBGMSG("TDG: solution size %d\n", feo->dh()->n_global_dofs());
273 
274  // Resize coefficient arrays
275  int qsize = max(feo->q<0>()->size(), max(feo->q<1>()->size(), max(feo->q<2>()->size(), feo->q<3>()->size())));
276  int max_edg_sides = max(mesh_->max_edge_sides(1), max(mesh_->max_edge_sides(2), mesh_->max_edge_sides(3)));
277  mm_coef.resize(qsize);
278  ad_coef.resize(n_subst_);
279  dif_coef.resize(n_subst_);
280  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
281  {
282  ad_coef[sbi].resize(qsize);
283  dif_coef[sbi].resize(qsize);
284  }
285  ad_coef_edg.resize(max_edg_sides);
286  dif_coef_edg.resize(max_edg_sides);
287  for (int sd=0; sd<max_edg_sides; sd++)
288  {
289  ad_coef_edg[sd].resize(n_subst_);
290  dif_coef_edg[sd].resize(n_subst_);
291  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
292  {
293  ad_coef_edg[sd][sbi].resize(qsize);
294  dif_coef_edg[sd][sbi].resize(qsize);
295  }
296  }
297 
298  // register output fields
299  output_rec = in_rec.val<Input::Record>("output_stream");
300  output_vec.resize(n_subst_);
301  output_solution.resize(n_subst_);
302  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
303  {
304  // for each substance we allocate output array and vector
305  output_solution[sbi] = new double[feo->dh()->n_global_dofs()];
306  VecCreateSeqWithArray(PETSC_COMM_SELF, 1, feo->dh()->n_global_dofs(), output_solution[sbi], &output_vec[sbi]);
307  }
308  data_.output_field.init(subst_names_);
309  data_.output_field.set_mesh(*mesh_);
310  data_.output_type(OutputTime::CORNER_DATA);
311 
312  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
313  {
314  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
315  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar> > output_field_ptr(new FieldFE<3, FieldValue<3>::Scalar>);
316  output_field_ptr->set_fe_data(feo->dh(), feo->mapping<1>(), feo->mapping<2>(), feo->mapping<3>(), &output_vec[sbi]);
317  data_.output_field[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
318  }
319  data_.set_limit_side(LimitSide::left);
322 
323  // set time marks for writing the output
325 
326  // allocate matrix and vector structures
327  ls = new LinSys*[n_subst_];
328  ls_dt = new LinSys_PETSC(feo->dh()->distr());
329  ( (LinSys_PETSC *)ls_dt )->set_from_input( in_rec.val<Input::Record>("solver") );
330  for (unsigned int sbi = 0; sbi < n_subst_; sbi++) {
331  ls[sbi] = new LinSys_PETSC(feo->dh()->distr());
332  ( (LinSys_PETSC *)ls[sbi] )->set_from_input( in_rec.val<Input::Record>("solver") );
333  ls[sbi]->set_solution(NULL);
334  }
335  stiffness_matrix = new Mat[n_subst_];
336  rhs = new Vec[n_subst_];
337 
338 }
339 
340 template<class Model>
342 {
343  delete time_;
344  delete ls_dt;
345 
346  if (feo->dh()->el_ds()->myp() == 0)
347  {
348  for (unsigned int i=0; i<n_subst_; i++)
349  {
350  VecDestroy(&output_vec[i]);
351  delete[] output_solution[i];
352  }
353  }
354 
355  for (unsigned int i=0; i<n_subst_; i++)
356  {
357  delete ls[i];
358  MatDestroy(&stiffness_matrix[i]);
359  VecDestroy(&rhs[i]);
360  }
361  delete[] ls;
362  delete[] stiffness_matrix;
363  delete[] rhs;
364  delete feo;
365  if (mass_balance_ != NULL) delete mass_balance_;
366 
367  gamma.clear();
368  subst_names_.clear();
369  delete output_stream;
370 }
371 
372 
373 template<class Model>
375 {
376  IS is;
377  VecScatter output_scatter;
378  int idx[] = { 0 };
379  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
380  {
381  // gather solution to output_vec[sbi]
382  ISCreateBlock(PETSC_COMM_SELF, ls[sbi]->size(), 1, idx, PETSC_COPY_VALUES, &is);
383  VecScatterCreate(ls[sbi]->get_solution(), is, output_vec[sbi], PETSC_NULL, &output_scatter);
384  VecScatterBegin(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
385  VecScatterEnd(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
386  VecScatterDestroy(&(output_scatter));
387  ISDestroy(&(is));
388  }
389 }
390 
391 
392 
393 template<class Model>
395 {
396  data_.set_time(*time_);
397 
398  // set initial conditions
399  set_initial_condition();
400  for (unsigned int sbi = 0; sbi < n_subst_; sbi++)
401  ( (LinSys_PETSC *)ls[sbi] )->set_initial_guess_nonzero();
402 
403  output_data();
404 }
405 
406 
407 template<class Model>
409 {
410  START_TIMER("DG-ONE STEP");
411 
412  time_->next_time();
413  time_->view("TDG");
414 
415  START_TIMER("data reinit");
416  data_.set_time(*time_);
417  END_TIMER("data reinit");
418 
419  // check first time assembly - needs preallocation
420  if (!allocation_done)
421  {
422  // preallocate mass matrix
423  ls_dt->start_allocation();
424  assemble_mass_matrix();
425  mass_matrix = NULL;
426 
427  // preallocate system matrix
428  for (unsigned int i=0; i<n_subst_; i++)
429  {
430  ls[i]->start_allocation();
431  stiffness_matrix[i] = NULL;
432  rhs[i] = NULL;
433  }
434  assemble_stiffness_matrix();
435  set_sources();
436  set_boundary_conditions();
437 
438  allocation_done = true;
439  }
440 
441  // assemble mass matrix
442  if (mass_matrix == NULL || data_.subset(FieldFlag::in_time_term).changed() )
443  {
444  ls_dt->start_add_assembly();
445  ls_dt->mat_zero_entries();
446  assemble_mass_matrix();
447  ls_dt->finish_assembly();
448  mass_matrix = *(ls_dt->get_matrix());
449  }
450 
451  // assemble stiffness matrix
452  if (stiffness_matrix[0] == NULL
453  || data_.subset(FieldFlag::in_main_matrix).changed()
454  || Model::flux_changed)
455  {
456  // new fluxes can change the location of Neumann boundary,
457  // thus stiffness matrix must be reassembled
458  for (unsigned int i=0; i<n_subst_; i++)
459  {
460  ls[i]->start_add_assembly();
461  ls[i]->mat_zero_entries();
462  }
463  assemble_stiffness_matrix();
464  for (unsigned int i=0; i<n_subst_; i++)
465  {
466  ls[i]->finish_assembly();
467 
468  if (stiffness_matrix[i] == NULL)
469  MatConvert(*( ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
470  else
471  MatCopy(*( ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
472  }
473  }
474 
475  // assemble right hand side (due to sources and boundary conditions)
476  if (rhs[0] == NULL
477  || data_.subset(FieldFlag::in_rhs).changed()
478  || Model::flux_changed)
479  {
480  for (unsigned int i=0; i<n_subst_; i++)
481  {
482  ls[i]->start_add_assembly();
483  ls[i]->rhs_zero_entries();
484  }
485  set_sources();
486  set_boundary_conditions();
487  for (unsigned int i=0; i<n_subst_; i++)
488  {
489  ls[i]->finish_assembly();
490 
491  VecDuplicate(*( ls[i]->get_rhs() ), &rhs[i]);
492  VecCopy(*( ls[i]->get_rhs() ), rhs[i]);
493  }
494  }
495 
496  Model::flux_changed = false;
497 
498 
499  /* Apply backward Euler time integration.
500  *
501  * Denoting A the stiffness matrix and M the mass matrix, the algebraic system at the k-th time level reads
502  *
503  * (1/dt M + A)u^k = f + 1/dt M.u^{k-1}
504  *
505  * Hence we modify at each time level the right hand side:
506  *
507  * f^k = f + 1/dt M u^{k-1},
508  *
509  * where f stands for the term stemming from the force and boundary conditions.
510  * Accordingly, we set
511  *
512  * A^k = A + 1/dt M.
513  *
514  */
515  Mat m;
516  START_TIMER("solve");
517  for (unsigned int i=0; i<n_subst_; i++)
518  {
519  MatConvert(stiffness_matrix[i], MATSAME, MAT_INITIAL_MATRIX, &m);
520  MatAXPY(m, 1./time_->dt(), mass_matrix, SUBSET_NONZERO_PATTERN);
521  ls[i]->set_matrix(m, DIFFERENT_NONZERO_PATTERN);
522  Vec y,w;
523  VecDuplicate(rhs[i], &y);
524  VecDuplicate(rhs[i], &w);
525  MatMult(mass_matrix, ls[i]->get_solution(), y);
526  VecWAXPY(w, 1./time_->dt(), y, rhs[i]);
527  ls[i]->set_rhs(w);
528 
529  VecDestroy(&y);
530  VecDestroy(&w);
531  MatDestroy(&m);
532 
533  ls[i]->solve();
534  }
535  END_TIMER("solve");
536 
537  if (mass_balance() != NULL)
538  mass_balance()->calculate(time_->t());
539 
540  END_TIMER("DG-ONE STEP");
541 }
542 
543 
544 template<class Model>
546 {
547  // So far the velocity_vector contains zeros, so we ignore it.
548  // Instead we use the value Side.flux.
549 
550  mh_dh = &dh;
551  Model::flux_changed = true;
552 
553 }
554 
555 
556 template<class Model>
558 {
559  if (!time_->is_current( time_->marks().type_output() )) return;
560 
561  START_TIMER("DG-OUTPUT");
562 
563  // gather the solution from all processors
564  output_vector_gather();
565  data_.subset(FieldFlag::allow_output).set_time( *time_);
566  data_.output(output_stream);
567  output_stream->write_time_frame();
568 
569  if (mass_balance() != NULL)
570  mass_balance()->output(time_->t());
571 
572  END_TIMER("DG-OUTPUT");
573 }
574 
575 
576 template<class Model>
578 {
579  START_TIMER("assemble_mass");
580  assemble_mass_matrix<1>();
581  assemble_mass_matrix<2>();
582  assemble_mass_matrix<3>();
583  END_TIMER("assemble_mass");
584 }
585 
586 
587 template<class Model> template<unsigned int dim>
589 {
590  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(), update_values | update_JxW_values | update_quadrature_points);
591  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
592  unsigned int dof_indices[ndofs];
593  PetscScalar local_mass_matrix[ndofs*ndofs];
594 
595  // assemble integral over elements
596  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
597  {
598  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
599  if (cell->dim() != dim) continue;
600 
601  fe_values.reinit(cell);
602  feo->dh()->get_dof_indices(cell, dof_indices);
603  ElementAccessor<3> ele_acc = cell->element_accessor();
604 
605  Model::compute_mass_matrix_coefficient(fe_values.point_list(), ele_acc, mm_coef);
606 
607  // assemble the local mass matrix
608  for (unsigned int i=0; i<ndofs; i++)
609  {
610  for (unsigned int j=0; j<ndofs; j++)
611  {
612  local_mass_matrix[i*ndofs+j] = 0;
613  for (unsigned int k=0; k<qsize; k++)
614  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);
615  }
616  }
617 
618  ls_dt->mat_set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, local_mass_matrix);
619  }
620 }
621 
622 
623 
624 
625 template<class Model>
627 {
628  START_TIMER("assemble_stiffness");
629  START_TIMER("assemble_volume_integrals");
630  assemble_volume_integrals<1>();
631  assemble_volume_integrals<2>();
632  assemble_volume_integrals<3>();
633  END_TIMER("assemble_volume_integrals");
634 
635  START_TIMER("assemble_fluxes_boundary");
636  assemble_fluxes_boundary<1>();
637  assemble_fluxes_boundary<2>();
638  assemble_fluxes_boundary<3>();
639  END_TIMER("assemble_fluxes_boundary");
640 
641  START_TIMER("assemble_fluxes_elem_elem");
642  assemble_fluxes_element_element<1>();
643  assemble_fluxes_element_element<2>();
644  assemble_fluxes_element_element<3>();
645  END_TIMER("assemble_fluxes_elem_elem");
646 
647  START_TIMER("assemble_fluxes_elem_side");
648  assemble_fluxes_element_side<1>();
649  assemble_fluxes_element_side<2>();
650  assemble_fluxes_element_side<3>();
651  END_TIMER("assemble_fluxes_elem_side");
652  END_TIMER("assemble_stiffness");
653 }
654 
655 
656 
657 template<class Model>
658 template<unsigned int dim>
660 {
661  FEValues<dim,3> fv_rt(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe_rt<dim>(),
663  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
665  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
666  unsigned int dof_indices[ndofs];
667  vector<arma::vec3> velocity(qsize);
668  vector<arma::vec> sources_sigma(qsize, arma::vec(n_substances()));
669  PetscScalar local_matrix[ndofs*ndofs];
670 
671  // assemble integral over elements
672  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
673  {
674  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
675  if (cell->dim() != dim) continue;
676 
677  fe_values.reinit(cell);
678  fv_rt.reinit(cell);
679  ElementAccessor<3> ele_acc = cell->element_accessor();
680  feo->dh()->get_dof_indices(cell, dof_indices);
681 
682  calculate_velocity(cell, velocity, fv_rt);
683  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), velocity, ele_acc, ad_coef, dif_coef);
684  Model::compute_sources_sigma(fe_values.point_list(), ele_acc, sources_sigma);
685 
686  // assemble the local stiffness matrix
687  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
688  {
689  for (unsigned int i=0; i<ndofs; i++)
690  for (unsigned int j=0; j<ndofs; j++)
691  local_matrix[i*ndofs+j] = 0;
692 
693  for (unsigned int k=0; k<qsize; k++)
694  {
695  for (unsigned int i=0; i<ndofs; i++)
696  {
697  arma::vec3 Kt_grad_i = dif_coef[sbi][k].t()*fe_values.shape_grad(i,k);
698  double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.shape_grad(i,k));
699 
700  for (unsigned int j=0; j<ndofs; j++)
701  local_matrix[i*ndofs+j] += (arma::dot(Kt_grad_i, fe_values.shape_grad(j,k))
702  -fe_values.shape_value(j,k)*ad_dot_grad_i
703  +sources_sigma[k][sbi]*fe_values.shape_value(j,k)*fe_values.shape_value(i,k))*fe_values.JxW(k);
704  }
705  }
706  ls[sbi]->mat_set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, local_matrix);
707  }
708  }
709 }
710 
711 
712 template<class Model>
714 {
715  START_TIMER("assemble_sources");
716  set_sources<1>();
717  set_sources<2>();
718  set_sources<3>();
719  END_TIMER("assemble_sources");
720 }
721 
722 template<class Model>
723 template<unsigned int dim>
725 {
726  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
728  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
729  vector<arma::vec> sources_conc(qsize, arma::vec(n_substances())),
730  sources_density(qsize, arma::vec(n_substances())),
731  sources_sigma(qsize, arma::vec(n_substances()));
732  unsigned int dof_indices[ndofs];
733  PetscScalar local_rhs[ndofs];
734  double source;
735 
736  // assemble integral over elements
737  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
738  {
739  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
740  if (cell->dim() != dim) continue;
741 
742  fe_values.reinit(cell);
743  feo->dh()->get_dof_indices(cell, dof_indices);
744 
745  Model::compute_source_coefficients(fe_values.point_list(), cell->element_accessor(), sources_conc, sources_density, sources_sigma);
746 
747  // assemble the local stiffness matrix
748  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
749  {
750  for (unsigned int i=0; i<ndofs; i++)
751  local_rhs[i] = 0;
752 
753  // compute sources
754  for (unsigned int k=0; k<qsize; k++)
755  {
756  source = (sources_density[k][sbi] + sources_conc[k][sbi]*sources_sigma[k][sbi])*fe_values.JxW(k);
757 
758  for (unsigned int i=0; i<ndofs; i++)
759  local_rhs[i] += source*fe_values.shape_value(i,k);
760  }
761  ls[sbi]->rhs_set_values(ndofs, (int *)dof_indices, local_rhs);
762  }
763  }
764 }
765 
766 
767 
768 template<class Model>
769 template<unsigned int dim>
771 {
772  vector<FESideValues<dim,3>*> fe_values;
773  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
774  update_values);
775  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size(),
776  n_max_sides = ad_coef_edg.size();
777  vector<unsigned int*> side_dof_indices;
778  PetscScalar local_matrix[ndofs*ndofs];
779  vector<vector<arma::vec3> > side_velocity(n_max_sides);
780  vector<arma::vec> dg_penalty(n_max_sides);
781  double gamma_l, omega[2], transport_flux;
782 
783  for (unsigned int sid=0; sid<n_max_sides; sid++)
784  {
785  side_dof_indices.push_back(new unsigned int[ndofs]);
786  fe_values.push_back(new FESideValues<dim,3>(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
788  }
789 
790  // assemble integral over sides
791  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
792  {
793  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
794  if (edg->n_sides < 2 || edg->side(0)->element()->dim() != dim) continue;
795 
796  for (int sid=0; sid<edg->n_sides; sid++)
797  {
798  typename DOFHandlerBase::CellIterator cell = edg->side(sid)->element();
799  ElementAccessor<3> ele_acc = cell->element_accessor();
800  feo->dh()->get_dof_indices(cell, side_dof_indices[sid]);
801  fe_values[sid]->reinit(cell, edg->side(sid)->el_idx());
802  fsv_rt.reinit(cell, edg->side(sid)->el_idx());
803  calculate_velocity(cell, side_velocity[sid], fsv_rt);
804  Model::compute_advection_diffusion_coefficients(fe_values[sid]->point_list(), side_velocity[sid], ele_acc, ad_coef_edg[sid], dif_coef_edg[sid]);
805  dg_penalty[sid] = data_.dg_penalty.value(cell->centre(), ele_acc);
806  }
807 
808  // fluxes and penalty
809  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
810  {
811  vector<double> fluxes(edg->n_sides);
812  for (int sid=0; sid<edg->n_sides; sid++)
813  {
814  fluxes[sid] = 0;
815  for (unsigned int k=0; k<qsize; k++)
816  fluxes[sid] += arma::dot(ad_coef_edg[sid][sbi][k], fe_values[sid]->normal_vector(k))*fe_values[sid]->JxW(k);
817  fluxes[sid] /= edg->side(sid)->measure();
818  }
819 
820  for (int s1=0; s1<edg->n_sides; s1++)
821  {
822  for (int s2=s1+1; s2<edg->n_sides; s2++)
823  {
824  ASSERT(edg->side(s1)->valid(), "Invalid side of edge.");
825  if (!feo->dh()->el_is_local(edg->side(s1)->element().index())
826  && !feo->dh()->el_is_local(edg->side(s2)->element().index())) continue;
827 
828  arma::vec3 nv = fe_values[s1]->normal_vector(0);
829 
830  // set up the parameters for DG method
831  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);
832 
833  int sd[2];
834  sd[0] = s1;
835  sd[1] = s2;
836 
837 #define AVERAGE(i,k,side_id) (fe_values[sd[side_id]]->shape_value(i,k)*0.5)
838 #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])
839 #define JUMP(i,k,side_id) ((side_id==0?1:-1)*fe_values[sd[side_id]]->shape_value(i,k))
840 
841  // For selected pair of elements:
842  for (int n=0; n<2; n++)
843  {
844  if (!feo->dh()->el_is_local(edg->side(sd[n])->element().index())) continue;
845 
846  for (int m=0; m<2; m++)
847  {
848  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
849  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
850  local_matrix[i*fe_values[sd[m]]->n_dofs()+j] = 0;
851 
852  for (unsigned int k=0; k<qsize; k++)
853  {
854  double flux_times_JxW = transport_flux*fe_values[0]->JxW(k);
855  double gamma_times_JxW = gamma_l*fe_values[0]->JxW(k);
856 
857  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
858  {
859  double flux_JxW_jump_i = flux_times_JxW*JUMP(i,k,n);
860  double gamma_JxW_jump_i = gamma_times_JxW*JUMP(i,k,n);
861  double JxW_jump_i = fe_values[0]->JxW(k)*JUMP(i,k,n);
862  double JxW_var_wavg_i = fe_values[0]->JxW(k)*WAVERAGE(i,k,n)*dg_variant;
863 
864  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
865  {
866  int index = i*fe_values[sd[m]]->n_dofs()+j;
867 
868  // flux due to transport (applied on interior edges) (average times jump)
869  local_matrix[index] += flux_JxW_jump_i*AVERAGE(j,k,m);
870 
871  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
872  local_matrix[index] += gamma_JxW_jump_i*JUMP(j,k,m);
873 
874  // terms due to diffusion
875  local_matrix[index] -= WAVERAGE(j,k,m)*JxW_jump_i;
876  local_matrix[index] -= JUMP(j,k,m)*JxW_var_wavg_i;
877  }
878  }
879  }
880  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);
881  }
882  }
883 #undef AVERAGE
884 #undef WAVERAGE
885 #undef JUMP
886  }
887  }
888  }
889  }
890 
891  for (unsigned int i=0; i<n_max_sides; i++)
892  {
893  delete fe_values[i];
894  delete[] side_dof_indices[i];
895  }
896 }
897 
898 
899 template<class Model>
900 template<unsigned int dim>
902 {
903  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
905  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
906  update_values);
907  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
908  unsigned int side_dof_indices[ndofs];
909  PetscScalar local_matrix[ndofs*ndofs];
910  vector<arma::vec3> side_velocity;
911  vector<arma::vec> robin_sigma(qsize, arma::vec(n_substances()));
912  arma::vec dg_penalty;
913  double gamma_l;
914 
915  // assemble boundary integral
916  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
917  {
918  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
919  if (edg->n_sides > 1) continue;
920  // check spatial dimension
921  if (edg->side(0)->dim() != dim-1) continue;
922  // skip edges lying not on the boundary
923  if (edg->side(0)->cond() == NULL) continue;
924 
925  SideIter side = edg->side(0);
926  typename DOFHandlerBase::CellIterator cell = side->element();
927  ElementAccessor<3> ele_acc = cell->element_accessor();
928  feo->dh()->get_dof_indices(cell, side_dof_indices);
929  fe_values_side.reinit(cell, side->el_idx());
930  fsv_rt.reinit(cell, side->el_idx());
931 
932  calculate_velocity(cell, side_velocity, fsv_rt);
933  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), side_velocity, ele_acc, ad_coef, dif_coef);
934  dg_penalty = data_.dg_penalty.value(cell->centre(), ele_acc);
935  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), side->cond()->element_accessor());
936  data_.bc_robin_sigma.value_list(fe_values_side.point_list(), side->cond()->element_accessor(), robin_sigma);
937 
938  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
939  {
940  for (unsigned int i=0; i<ndofs; i++)
941  for (unsigned int j=0; j<ndofs; j++)
942  local_matrix[i*ndofs+j] = 0;
943 
944  // On Neumann boundaries we have only term from integrating by parts the advective term,
945  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
946  double side_flux = 0;
947  for (unsigned int k=0; k<qsize; k++)
948  side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
949  double transport_flux = side_flux/side->measure();
950 
951  if ((bc_type[sbi] == EqData::dirichlet) || (bc_type[sbi] == EqData::inflow ))
952  {
953  // set up the parameters for DG method
954  set_DG_parameters_boundary(side, qsize, dif_coef[sbi], transport_flux, fe_values_side.normal_vector(0), dg_penalty[sbi], gamma_l);
955  gamma[sbi][side->cond_idx()] = gamma_l;
956  if (bc_type[sbi] == EqData::dirichlet || side_flux < -mh_dh->precision())
957  transport_flux += gamma_l;
958  }
959 
960  // fluxes and penalty
961  for (unsigned int k=0; k<qsize; k++)
962  {
963  double flux_times_JxW;
964  if (bc_type[sbi] == EqData::robin)
965  flux_times_JxW = (transport_flux + robin_sigma[k][sbi])*fe_values_side.JxW(k);
966  else
967  flux_times_JxW = transport_flux*fe_values_side.JxW(k);
968 
969  for (unsigned int i=0; i<ndofs; i++)
970  for (unsigned int j=0; j<ndofs; j++)
971  {
972  // flux due to advection and penalty
973  local_matrix[i*ndofs+j] += flux_times_JxW*fe_values_side.shape_value(i,k)*fe_values_side.shape_value(j,k);
974 
975  // flux due to diffusion (only on dirichlet and inflow boundary)
976  if (bc_type[sbi] == EqData::dirichlet || (bc_type[sbi] == EqData::inflow && side_flux < -mh_dh->precision()))
977  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)
978  + 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
979  )*fe_values_side.JxW(k);
980  }
981  }
982  ls[sbi]->mat_set_values(ndofs, (int *)side_dof_indices, ndofs, (int *)side_dof_indices, local_matrix);
983  }
984  }
985 }
986 
987 
988 template<class Model>
989 template<unsigned int dim>
991 {
992 
993  if (dim == 1) return;
994  FEValues<dim-1,3> fe_values_vb(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe<dim-1>(),
996  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
998  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
999  update_values);
1000  FEValues<dim-1,3> fv_rt(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe_rt<dim-1>(),
1001  update_values);
1002 
1003  vector<FEValuesSpaceBase<3>*> fv_sb(2);
1004  const unsigned int ndofs = feo->fe<dim>()->n_dofs(); // number of local dofs
1005  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
1006  unsigned int side_dof_indices[2*ndofs], n_dofs[2];
1007  vector<arma::vec3> velocity_higher, velocity_lower;
1008  vector<arma::vec> frac_sigma(qsize, arma::vec(n_substances()));
1009  vector<double> csection_lower(qsize), csection_higher(qsize), mm_coef_lower(qsize), mm_coef_higher(qsize);
1010  PetscScalar local_matrix[4*ndofs*ndofs];
1011  double comm_flux[2][2];
1012 
1013  // index 0 = element with lower dimension,
1014  // index 1 = side of element with higher dimension
1015  fv_sb[0] = &fe_values_vb;
1016  fv_sb[1] = &fe_values_side;
1017 
1018  // assemble integral over sides
1019  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
1020  {
1021  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
1022  // skip neighbours of different dimension
1023  if (nb->element()->dim() != dim-1) continue;
1024 
1025  typename DOFHandlerBase::CellIterator cell_sub = mesh_->element.full_iter(nb->element());
1026  feo->dh()->get_dof_indices(cell_sub, side_dof_indices);
1027  fe_values_vb.reinit(cell_sub);
1028  n_dofs[0] = fv_sb[0]->n_dofs();
1029 
1030  typename DOFHandlerBase::CellIterator cell = nb->side()->element();
1031  feo->dh()->get_dof_indices(cell, side_dof_indices+n_dofs[0]);
1032  fe_values_side.reinit(cell, nb->side()->el_idx());
1033  n_dofs[1] = fv_sb[1]->n_dofs();
1034 
1035  // Element id's for testing if they belong to local partition.
1036  int element_id[2];
1037  element_id[0] = cell_sub.index();
1038  element_id[1] = cell.index();
1039 
1040  fsv_rt.reinit(cell, nb->side()->el_idx());
1041  fv_rt.reinit(cell_sub);
1042  calculate_velocity(cell, velocity_higher, fsv_rt);
1043  calculate_velocity(cell_sub, velocity_lower, fv_rt);
1044  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_lower, cell_sub->element_accessor(), ad_coef_edg[0], dif_coef_edg[0]);
1045  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_higher, cell->element_accessor(), ad_coef_edg[1], dif_coef_edg[1]);
1046  Model::compute_mass_matrix_coefficient(fe_values_vb.point_list(), cell_sub->element_accessor(), mm_coef_lower);
1047  Model::compute_mass_matrix_coefficient(fe_values_vb.point_list(), cell->element_accessor(), mm_coef_higher);
1048  data_.cross_section.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), csection_lower);
1049  data_.cross_section.value_list(fe_values_vb.point_list(), cell->element_accessor(), csection_higher);
1050  data_.fracture_sigma.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), frac_sigma);
1051 
1052  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1053  {
1054  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
1055  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
1056  local_matrix[i*(n_dofs[0]+n_dofs[1])+j] = 0;
1057 
1058  // set transmission conditions
1059  for (unsigned int k=0; k<qsize; k++)
1060  {
1061  /* The communication flux has two parts:
1062  * - "diffusive" term containing sigma
1063  * - "advective" term representing usual upwind
1064  *
1065  * The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
1066  * than b and A in the manual.
1067  * In calculation of sigma there appears one more csection_lower in the denominator.
1068  */
1069  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))*
1070  2*csection_higher[k]*csection_higher[k]/(csection_lower[k]*csection_lower[k]);
1071 
1072  // Since mm_coef_* contains cross section, we have to divide by it.
1073  double transport_flux = arma::dot(ad_coef_edg[1][sbi][k], fe_values_side.normal_vector(k));
1074  double por_lower_over_higher = mm_coef_lower[k]*csection_higher[k]/(mm_coef_higher[k]*csection_lower[k]);
1075 
1076  comm_flux[0][0] = (sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1077  comm_flux[0][1] = -(sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1078  comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1079  comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1080 
1081  for (int n=0; n<2; n++)
1082  {
1083  if (!feo->dh()->el_is_local(element_id[n])) continue;
1084 
1085  for (unsigned int i=0; i<n_dofs[n]; i++)
1086  for (int m=0; m<2; m++)
1087  for (unsigned int j=0; j<n_dofs[m]; j++)
1088  local_matrix[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
1089  comm_flux[m][n]*fv_sb[m]->shape_value(j,k)*fv_sb[n]->shape_value(i,k);
1090  }
1091  }
1092  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);
1093  }
1094  }
1095 
1096 }
1097 
1098 
1099 
1100 
1101 
1102 
1103 template<class Model>
1105 {
1106  START_TIMER("assemble_bc");
1107  set_boundary_conditions<1>();
1108  set_boundary_conditions<2>();
1109  set_boundary_conditions<3>();
1110  END_TIMER("assemble_bc");
1111 }
1112 
1113 
1114 template<class Model>
1115 template<unsigned int dim>
1117 {
1118  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1120  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1121  update_values);
1122  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1123  unsigned int side_dof_indices[ndofs];
1124  double local_rhs[ndofs];
1125  vector<arma::vec> bc_values(qsize, arma::vec(n_substances())),
1126  bc_fluxes(qsize, arma::vec(n_substances())),
1127  bc_sigma(qsize, arma::vec(n_substances()));
1128  vector<arma::vec3> velocity;
1129  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1130  {
1131  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
1132  if (edg->n_sides > 1) continue;
1133  if (edg->side(0)->dim() != dim-1) continue;
1134  // skip edges lying not on the boundary
1135  if (edg->side(0)->cond() == NULL) continue;
1136 
1137  SideIter side = edg->side(0);
1138  typename DOFHandlerBase::CellIterator cell = mesh().element.full_iter(side->element());
1139  ElementAccessor<3> ele_acc = side->cond()->element_accessor();
1140 
1141  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), ele_acc);
1142 
1143  fe_values_side.reinit(cell, side->el_idx());
1144  fsv_rt.reinit(cell, side->el_idx());
1145  calculate_velocity(cell, velocity, fsv_rt);
1146 
1147  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), velocity, side->element()->element_accessor(), ad_coef, dif_coef);
1148  Model::compute_dirichlet_bc(fe_values_side.point_list(), ele_acc, bc_values);
1149  data_.bc_flux.value_list(fe_values_side.point_list(), ele_acc, bc_fluxes);
1150  data_.bc_robin_sigma.value_list(fe_values_side.point_list(), ele_acc, bc_sigma);
1151 
1152  feo->dh()->get_dof_indices(cell, side_dof_indices);
1153 
1154  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1155  {
1156  for (unsigned int i=0; i<ndofs; i++) local_rhs[i] = 0;
1157 
1158  for (unsigned int k=0; k<qsize; k++)
1159  {
1160  double bc_term = 0;
1161  arma::vec3 bc_grad;
1162  bc_grad.zeros();
1163  if ((bc_type[sbi] == EqData::inflow && mh_dh->side_flux( *edg->side(0) ) < -mh_dh->precision())
1164  || (bc_type[sbi] == EqData::dirichlet))
1165  {
1166  bc_term = gamma[sbi][side->cond_idx()]*bc_values[k][sbi]*fe_values_side.JxW(k);
1167  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));
1168  }
1169  else if (bc_type[sbi] == EqData::neumann)
1170  {
1171  bc_term = bc_fluxes[k][sbi]*fe_values_side.JxW(k);
1172  }
1173  else if (bc_type[sbi] == EqData::robin)
1174  {
1175  bc_term = bc_sigma[k][sbi]*bc_values[k][sbi]*fe_values_side.JxW(k);
1176  }
1177 
1178  for (unsigned int i=0; i<ndofs; i++)
1179  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k)
1180  + arma::dot(bc_grad,fe_values_side.shape_grad(i,k));
1181  }
1182  ls[sbi]->rhs_set_values(ndofs, (int *)side_dof_indices, local_rhs);
1183  }
1184  }
1185 }
1186 
1187 
1188 
1189 // TODO: The detection of side number from SideIter
1190 // in TransportDG::calculate_velocity()
1191 // should be done with help of class RefElement. This however requires
1192 // that the MH_DofHandler uses the node/side ordering defined in
1193 // the respective RefElement.
1194 template<class Model>
1195 template<unsigned int dim>
1197 {
1198  std::map<const Node*, int> node_nums;
1199  for (unsigned int i=0; i<cell->n_nodes(); i++)
1200  node_nums[cell->node[i]] = i;
1201 
1202  velocity.resize(fv.n_points());
1203 
1204  for (unsigned int k=0; k<fv.n_points(); k++)
1205  {
1206  velocity[k].zeros();
1207  for (unsigned int sid=0; sid<cell->n_sides(); sid++)
1208  {
1209  if (cell->side(sid)->dim() != dim-1) continue;
1210  int num = dim*(dim+1)/2;
1211  for (unsigned int i=0; i<cell->side(sid)->n_nodes(); i++)
1212  num -= node_nums[cell->side(sid)->node(i)];
1213  velocity[k] += fv.shape_vector(num,k) * mh_dh->side_flux( *(cell->side(sid)) );
1214  }
1215  }
1216 }
1217 
1218 
1219 
1220 
1221 
1222 template<class Model>
1224  const int s1,
1225  const int s2,
1226  const int K_size,
1227  const vector<arma::mat33> &K1,
1228  const vector<arma::mat33> &K2,
1229  const vector<double> &fluxes,
1230  const arma::vec3 &normal_vector,
1231  const double alpha1,
1232  const double alpha2,
1233  double &gamma,
1234  double *omega,
1235  double &transport_flux)
1236 {
1237  double delta[2];
1238  double h = 0;
1239  double local_alpha = 0;
1240 
1241  ASSERT(edg.side(s1)->valid(), "Invalid side of an edge.");
1242  SideIter s = edg.side(s1);
1243 
1244  // calculate the side diameter
1245  if (s->dim() == 0)
1246  {
1247  h = 1;
1248  }
1249  else
1250  {
1251  for (unsigned int i=0; i<s->n_nodes(); i++)
1252  for (unsigned int j=i+1; j<s->n_nodes(); j++)
1253  h = max(h, s->node(i)->distance(*s->node(j)));
1254  }
1255 
1256  // calculate the total in- and out-flux through the edge
1257  double pflux = 0, nflux = 0;
1258  for (int i=0; i<edg.n_sides; i++)
1259  {
1260  if (fluxes[i] > 0)
1261  pflux += fluxes[i];
1262  else
1263  nflux += fluxes[i];
1264  }
1265 
1266  // calculate the flux from s1 to s2
1267  if (fluxes[s2] > 0 && fluxes[s1] < 0 && s1 < s2)
1268  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
1269  else if (fluxes[s2] < 0 && fluxes[s1] > 0 && s1 < s2)
1270  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
1271  else if (s1==s2)
1272  transport_flux = fluxes[s1];
1273  else
1274  transport_flux = 0;
1275 
1276  gamma = 0.5*fabs(transport_flux);
1277 
1278 
1279  // determine local DG penalty
1280  local_alpha = max(alpha1, alpha2);
1281 
1282  if (s1 == s2)
1283  {
1284  omega[0] = 1;
1285 
1286  // delta is set to the average value of Kn.n on the side
1287  delta[0] = 0;
1288  for (int k=0; k<K_size; k++)
1289  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1290  delta[0] /= K_size;
1291 
1292  gamma += local_alpha/h*(delta[0]);
1293  }
1294  else
1295  {
1296  delta[0] = 0;
1297  delta[1] = 0;
1298  for (int k=0; k<K_size; k++)
1299  {
1300  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1301  delta[1] += dot(K2[k]*normal_vector,normal_vector);
1302  }
1303  delta[0] /= K_size;
1304  delta[1] /= K_size;
1305 
1306  double delta_sum = delta[0] + delta[1];
1307 
1308  if (delta_sum > numeric_limits<double>::epsilon())
1309  {
1310  omega[0] = delta[1]/delta_sum;
1311  omega[1] = delta[0]/delta_sum;
1312  gamma += local_alpha/h*(delta[0]*delta[1]/delta_sum);
1313  }
1314  else
1315  for (int i=0; i<2; i++) omega[i] = 0;
1316  }
1317 }
1318 
1319 
1320 
1321 
1322 
1323 
1324 template<class Model>
1326  const int K_size,
1327  const vector<arma::mat33> &K,
1328  const double flux,
1329  const arma::vec3 &normal_vector,
1330  const double alpha,
1331  double &gamma)
1332 {
1333  double delta = 0, h = 0;
1334 
1335  // calculate the side diameter
1336  if (side->dim() == 0)
1337  {
1338  h = 1;
1339  }
1340  else
1341  {
1342  for (unsigned int i=0; i<side->n_nodes(); i++)
1343  for (unsigned int j=i+1; j<side->n_nodes(); j++)
1344  h = max(h, side->node(i)->distance( *side->node(j) ));
1345  }
1346 
1347  // delta is set to the average value of Kn.n on the side
1348  for (int k=0; k<K_size; k++)
1349  delta += dot(K[k]*normal_vector,normal_vector);
1350  delta /= K_size;
1351 
1352  gamma = 0.5*fabs(flux) + alpha/h*delta;
1353 }
1354 
1355 
1356 
1357 
1358 
1359 template<class Model>
1361 {
1362  START_TIMER("set_init_cond");
1363  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1364  ls[sbi]->start_allocation();
1365  prepare_initial_condition<1>();
1366  prepare_initial_condition<2>();
1367  prepare_initial_condition<3>();
1368 
1369  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1370  ls[sbi]->start_add_assembly();
1371  prepare_initial_condition<1>();
1372  prepare_initial_condition<2>();
1373  prepare_initial_condition<3>();
1374 
1375  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1376  {
1377  ls[sbi]->finish_assembly();
1378  ls[sbi]->solve();
1379  }
1380  END_TIMER("set_init_cond");
1381 }
1382 
1383 template<class Model>
1384 template<unsigned int dim>
1386 {
1387  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1389  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1390  unsigned int dof_indices[ndofs];
1391  double matrix[ndofs*ndofs], rhs[ndofs];
1392  std::vector<arma::vec> init_values(qsize);
1393 
1394  for (unsigned int k=0; k<qsize; k++)
1395  init_values[k].resize(n_subst_);
1396 
1397  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
1398  {
1399  typename DOFHandlerBase::CellIterator elem = mesh_->element(feo->dh()->el_index(i_cell));
1400  if (elem->dim() != dim) continue;
1401 
1402  ElementAccessor<3> ele_acc = elem->element_accessor();
1403  feo->dh()->get_dof_indices(elem, dof_indices);
1404  fe_values.reinit(elem);
1405 
1406  Model::compute_init_cond(fe_values.point_list(), ele_acc, init_values);
1407 
1408  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1409  {
1410  for (unsigned int i=0; i<ndofs; i++)
1411  {
1412  rhs[i] = 0;
1413  for (unsigned int j=0; j<ndofs; j++)
1414  matrix[i*ndofs+j] = 0;
1415  }
1416 
1417  for (unsigned int k=0; k<qsize; k++)
1418  {
1419  double rhs_term = init_values[k](sbi)*fe_values.JxW(k);
1420 
1421  for (unsigned int i=0; i<ndofs; i++)
1422  {
1423  for (unsigned int j=0; j<ndofs; j++)
1424  matrix[i*ndofs+j] += fe_values.shape_value(i,k)*fe_values.shape_value(j,k)*fe_values.JxW(k);
1425 
1426  rhs[i] += fe_values.shape_value(i,k)*rhs_term;
1427  }
1428  }
1429  ls[sbi]->set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, matrix, rhs);
1430  }
1431  }
1432 }
1433 
1434 
1435 template<class Model>
1436 void TransportDG<Model>::calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance)
1437 {
1438  calc_fluxes<1>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1439  calc_fluxes<2>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1440  calc_fluxes<3>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1441 }
1442 
1443 template<class Model>
1444 template<unsigned int dim>
1445 void TransportDG<Model>::calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance)
1446 {
1447  FESideValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1449  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1450  update_values);
1451  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1452  unsigned int dof_indices[ndofs];
1453  vector<arma::vec3> side_velocity(qsize);
1454  vector<arma::vec> bc_values(qsize);
1455  arma::vec3 conc_grad;
1456 
1457  for (unsigned int k=0; k<qsize; k++)
1458  bc_values[k].resize(n_subst_);
1459 
1460  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1461  {
1462  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
1463  if (edg->n_sides > 1) continue;
1464  if (edg->side(0)->dim() != dim-1) continue;
1465  // skip edges lying not on the boundary
1466  if (edg->side(0)->cond() == NULL) continue;
1467 
1468  SideIter side = edg->side(0);
1470  ElementAccessor<3> ele_acc = cell->element_accessor();
1471  Region r = side->cond()->element_accessor().region();
1472  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d, edg: % d\n", side->cond()->bc_ele_idx_, side->cond()->edge_idx_);
1473  unsigned int bc_region_idx = r.boundary_idx();
1474 
1475  fe_values.reinit(cell, side->el_idx());
1476  fsv_rt.reinit(cell, side->el_idx());
1477  feo->dh()->get_dof_indices(cell, dof_indices);
1478  calculate_velocity(cell, side_velocity, fsv_rt);
1479  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), side_velocity, ele_acc, ad_coef, dif_coef);
1480  Model::compute_dirichlet_bc(fe_values.point_list(), side->cond()->element_accessor(), bc_values);
1481  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), side->cond()->element_accessor());
1482 
1483  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1484  {
1485  double mass_flux = 0;
1486  double water_flux = 0;
1487  for (unsigned int k=0; k<qsize; k++)
1488  water_flux += arma::dot(ad_coef[sbi][k], fe_values.normal_vector(k))*fe_values.JxW(k);
1489  water_flux /= side->measure();
1490 
1491  for (unsigned int k=0; k<qsize; k++)
1492  {
1493  double conc = 0;
1494  conc_grad.zeros();
1495  for (unsigned int i=0; i<ndofs; i++)
1496  {
1497  conc += fe_values.shape_value(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1498  conc_grad += fe_values.shape_grad(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1499  }
1500 
1501  // flux due to advection
1502  mass_flux += water_flux*conc*fe_values.JxW(k);
1503 
1504  if (bc_type[sbi] == EqData::dirichlet || (bc_type[sbi] == EqData::inflow && water_flux*side->measure() < -mh_dh->precision()))
1505  {
1506  // flux due to diffusion
1507  mass_flux -= arma::dot(dif_coef[sbi][k]*conc_grad,fe_values.normal_vector(k))*fe_values.JxW(k);
1508 
1509  // the penalty term has to be added otherwise the mass balance will not hold
1510  mass_flux -= gamma[sbi][side->cond_idx()]*(bc_values[k][sbi] - conc)*fe_values.JxW(k);
1511  }
1512 
1513  }
1514  bcd_balance[sbi][bc_region_idx] += mass_flux;
1515  if (mass_flux >= 0) bcd_plus_balance[sbi][bc_region_idx] += mass_flux;
1516  else bcd_minus_balance[sbi][bc_region_idx] += mass_flux;
1517  }
1518  }
1519 
1520 }
1521 
1522 template<class Model>
1524 {
1525  calc_elem_sources<1>(mass, src_balance);
1526  calc_elem_sources<2>(mass, src_balance);
1527  calc_elem_sources<3>(mass, src_balance);
1528 }
1529 
1530 template<class Model>
1531 template<unsigned int dim>
1533 {
1534  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1536  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1537  unsigned int dof_indices[ndofs];
1538  vector<arma::vec> sources_conc(qsize, arma::vec(n_substances())),
1539  sources_density(qsize, arma::vec(n_substances())),
1540  sources_sigma(qsize, arma::vec(n_substances()));
1541  double mass_sum, sources_sum, conc, conc_diff;
1542 
1543  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
1544  {
1545  typename DOFHandlerBase::CellIterator elem = mesh_->element(feo->dh()->el_index(i_cell));
1546  if (elem->dim() != dim) continue;
1547 
1548  Region r = elem->element_accessor().region();
1549  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d\n", elem.index());
1550  unsigned int region_idx = r.bulk_idx();
1551 
1552  ElementAccessor<3> ele_acc = elem->element_accessor();
1553  fe_values.reinit(elem);
1554  feo->dh()->get_dof_indices(elem, dof_indices);
1555 
1556  Model::compute_mass_matrix_coefficient(fe_values.point_list(), ele_acc, mm_coef);
1557  Model::compute_source_coefficients(fe_values.point_list(), ele_acc, sources_conc, sources_density, sources_sigma);
1558 
1559  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1560  {
1561  mass_sum = 0;
1562  sources_sum = 0;
1563 
1564  for (unsigned int k=0; k<qsize; k++)
1565  {
1566  conc = 0;
1567  for (unsigned int i=0; i<ndofs; i++)
1568  conc += fe_values.shape_value(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1569 
1570  mass_sum += mm_coef[k]*conc*fe_values.JxW(k);
1571 
1572  conc_diff = sources_conc[k][sbi] - conc;
1573  sources_sum += (sources_density[k][sbi] + conc_diff*sources_sigma[k][sbi])*fe_values.JxW(k);
1574  }
1575 
1576  mass[sbi][region_idx] += mass_sum;
1577  src_balance[sbi][region_idx] += sources_sum;
1578  }
1579  }
1580 }
1581 
1582 
1583 
1584 
1585 
1587 template class TransportDG<HeatTransferModel>;
1588 
1589 
1590 
1591 
static Input::Type::Selection output_selection
void calculate(double time)
Definition: mass_balance.cc:89
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
#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.
vector< double > mm_coef
Mass matrix coefficients.
Header: The functions for all outputs.
Accessor to input data conforming to declared Array.
Definition: accessors.hh:521
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
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:39
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
double distance(const Node &n2) const
Definition: nodes.hh:98
void assemble_mass_matrix()
Assembles the mass matrix.
Boundary * cond() const
Definition: side_impl.hh:59
double fix_dt_until_mark()
Fixing time step until fixed time mark.
OutputTime * output_stream
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.
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
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
TimeMark::Type type_output()
Definition: time_marks.hh:192
static Default obligatory()
Definition: type_record.hh:87
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
MassBalance * mass_balance_
object for calculation and writing the mass balance to file.
Definition: mesh.h:108
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...
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:148
void set_from_input(const Input::Record in_rec)
#define AVERAGE(i, k, side_id)
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:341
double t() const
vector< Boundary > boundary_
Definition: mesh.h:202
Basic time management functionality for unsteady (and steady) solvers (class Equation).
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
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
static TimeMarks & marks()
Class for declaration of inputs sequences.
Definition: type_base.hh:230
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)
const unsigned int n_points()
Returns the number of quadrature points.
Definition: fe_values.cc:202
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
Definition: region.hh:75
double precision() const
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 write all registered data to output streams.
Definition: output.cc:171
unsigned int dim() const
Definition: side_impl.hh:26
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
FEObjects * feo
Finite element objects.
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:90
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:44
Definitions of basic Lagrangean finite elements with polynomial shape functions.
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:308
FieldCommon & input_selection(const Input::Type::Selection *element_selection)
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:100
Mat mass_matrix
The mass matrix.
#define START_TIMER(tag)
Starts a timer with specified tag.
FiniteElement< dim, spacedim > * fe
The used finite element.
Definition: fe_values.hh:337
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:251
Mat * stiffness_matrix
The stiffness matrix.
void output(double time)
Write computed fields to file.
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).
unsigned int dg_order
Polynomial order of finite elements.
void output_vector_gather()
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh: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.cc:331
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
DOFHandlerMultiDim * dh()
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 vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Definition: fe_values.cc:190
std::vector< string > subst_names_
Names of transported substances.
const double epsilon
Definition: mathfce.h:6
void output_data() override
Write computed fields.
vector< Neighbour > vb_neighbours_
Definition: mesh.h:228
double side_flux(const Side &side) const
temporary replacement for DofHandler accessor, flux through given side
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
FieldCommon & description(const string &description)
Definition: field_common.hh:78
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:175
Discontinuous Galerkin method for equation of transport with dispersion.
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
void add_admissible_field_names(const Input::Array &in_array, const Input::Type::Selection &in_sel)
Registers names of output fields that can be written using this stream.
Definition: output.cc:231
arma::vec3 centre() const
Definition: elements.cc:130
Definition: system.hh:72
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:312
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
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
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:209
FullIter full_iter(Iter it)
Definition: sys_vector.hh:438
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:71
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
unsigned int n_substances() override
Returns number of trnasported substances.
Record type proxy class.
Definition: type_record.hh:161
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
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:40
Base class for FEValues and FESideValues.
Definition: fe_values.hh:43
Class for representation SI units of Fields.
Definition: unit_si.hh:31
const Node * node(unsigned int i) const
Definition: side_impl.hh:35
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.
ElementAccessor< 3 > element_accessor()
Definition: boundaries.cc:55
friend class MassBalance
Definition: mass_balance.hh:62
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)
Definition: transport_dg.cc:84
TimeGovernor * time_
Definition: equation.hh:219
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:198
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