Flow123d  jenkins-Flow123d-windows32-release-multijob-28
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).")
193  .input_default("1.0")
195 
196  *this+=dg_penalty
197  .name("dg_penalty")
198  .description(
199  "Penalty parameter influencing the discontinuity of the solution (for each substance). "
200  "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
201  .input_default("1.0")
203 
204  *this+=bc_type
205  .name("bc_type")
206  .description(
207  "Boundary condition type, possible values: inflow, dirichlet, neumann, robin.")
208  .input_default("\"inflow\"")
211 
212 // std::vector<FieldEnum> list; list.push_back(neumann);
213  *this+=bc_flux
214  .name("bc_flux")
215  .description("Flux in Neumann boundary condition.")
216  .input_default("0.0")
218 // bc_flux.disable_where(bc_type, { dirichlet, inflow });
219  *this+=bc_robin_sigma
220  .name("bc_robin_sigma")
221  .description("Conductivity coefficient in Robin boundary condition.")
222  .input_default("0.0")
224 // bc_robin_sigma.disable_where(bc_type, {dirichlet, inflow, neumann});
225 
226  // add all input fields to the output list
227 
228 }
229 
230 template<class Model>
232  : TransportBase(init_mesh, in_rec),
233  mass_matrix(0),
234  allocation_done(false)
235 {
236  // Check that Model is derived from AdvectionDiffusionModel.
237  static_assert(std::is_base_of<AdvectionDiffusionModel, Model>::value, "");
238 
239  this->eq_data_ = &data_;
240 
241  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"));
243 
244  // Read names of transported substances.
245  Model::set_component_names(subst_names_, in_rec);
246  n_subst_ = subst_names_.size();
247 
248  Input::Iterator<Input::Record> it = in_rec.find<Input::Record>("mass_balance");
249  if (it)
250  mass_balance_ = new MassBalance(this, *it);
251 
252  // Set up physical parameters.
253  data_.set_mesh(init_mesh);
254  data_.set_n_components(n_subst_);
255  data_.set_input_list( in_rec.val<Input::Array>("input_fields") );
256  data_.set_limit_side(LimitSide::left);
257 
258 
259  // DG variant and order
260  dg_variant = in_rec.val<DGVariant>("dg_variant");
261  dg_order = in_rec.val<unsigned int>("dg_order");
262 
263  // DG stabilization parameters on boundary edges
264  gamma.resize(n_subst_);
265  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
266  gamma[sbi].resize(mesh_->boundary_.size());
267 
268  // create finite element structures and distribute DOFs
269  feo = new FEObjects(mesh_, dg_order);
270  DBGMSG("TDG: solution size %d\n", feo->dh()->n_global_dofs());
271 
272  // Resize coefficient arrays
273  int qsize = max(feo->q<0>()->size(), max(feo->q<1>()->size(), max(feo->q<2>()->size(), feo->q<3>()->size())));
274  int max_edg_sides = max(mesh_->max_edge_sides(1), max(mesh_->max_edge_sides(2), mesh_->max_edge_sides(3)));
275  mm_coef.resize(qsize);
276  ad_coef.resize(n_subst_);
277  dif_coef.resize(n_subst_);
278  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
279  {
280  ad_coef[sbi].resize(qsize);
281  dif_coef[sbi].resize(qsize);
282  }
283  ad_coef_edg.resize(max_edg_sides);
284  dif_coef_edg.resize(max_edg_sides);
285  for (int sd=0; sd<max_edg_sides; sd++)
286  {
287  ad_coef_edg[sd].resize(n_subst_);
288  dif_coef_edg[sd].resize(n_subst_);
289  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
290  {
291  ad_coef_edg[sd][sbi].resize(qsize);
292  dif_coef_edg[sd][sbi].resize(qsize);
293  }
294  }
295 
296  // register output fields
297  output_rec = in_rec.val<Input::Record>("output_stream");
298  output_vec.resize(n_subst_);
299  output_solution.resize(n_subst_);
300  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
301  {
302  // for each substance we allocate output array and vector
303  output_solution[sbi] = new double[feo->dh()->n_global_dofs()];
304  VecCreateSeqWithArray(PETSC_COMM_SELF, 1, feo->dh()->n_global_dofs(), output_solution[sbi], &output_vec[sbi]);
305  }
306  data_.output_field.init(subst_names_);
307  data_.output_field.set_mesh(*mesh_);
308  data_.output_type(OutputTime::CORNER_DATA);
309 
310  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
311  {
312  // create shared pointer to a FieldFE, pass FE data and push this FieldFE to output_field on all regions
313  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar> > output_field_ptr(new FieldFE<3, FieldValue<3>::Scalar>);
314  output_field_ptr->set_fe_data(feo->dh(), feo->mapping<1>(), feo->mapping<2>(), feo->mapping<3>(), &output_vec[sbi]);
315  data_.output_field[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
316  }
317  data_.set_limit_side(LimitSide::left);
320 
321  // set time marks for writing the output
323 
324  // allocate matrix and vector structures
325  ls = new LinSys*[n_subst_];
326  ls_dt = new LinSys_PETSC(feo->dh()->distr());
327  ( (LinSys_PETSC *)ls_dt )->set_from_input( in_rec.val<Input::Record>("solver") );
328  for (unsigned int sbi = 0; sbi < n_subst_; sbi++) {
329  ls[sbi] = new LinSys_PETSC(feo->dh()->distr());
330  ( (LinSys_PETSC *)ls[sbi] )->set_from_input( in_rec.val<Input::Record>("solver") );
331  ls[sbi]->set_solution(NULL);
332  }
333  stiffness_matrix = new Mat[n_subst_];
334  rhs = new Vec[n_subst_];
335 
336 }
337 
338 template<class Model>
340 {
341  delete time_;
342  delete ls_dt;
343 
344  if (feo->dh()->el_ds()->myp() == 0)
345  {
346  for (unsigned int i=0; i<n_subst_; i++)
347  {
348  VecDestroy(&output_vec[i]);
349  delete[] output_solution[i];
350  }
351  }
352 
353  for (unsigned int i=0; i<n_subst_; i++)
354  {
355  delete ls[i];
356  MatDestroy(&stiffness_matrix[i]);
357  VecDestroy(&rhs[i]);
358  }
359  delete[] ls;
360  delete[] stiffness_matrix;
361  delete[] rhs;
362  delete feo;
363  if (mass_balance_ != NULL) delete mass_balance_;
364 
365  gamma.clear();
366  subst_names_.clear();
367  delete output_stream;
368 }
369 
370 
371 template<class Model>
373 {
374  IS is;
375  VecScatter output_scatter;
376  int idx[] = { 0 };
377  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
378  {
379  // gather solution to output_vec[sbi]
380  ISCreateBlock(PETSC_COMM_SELF, ls[sbi]->size(), 1, idx, PETSC_COPY_VALUES, &is);
381  VecScatterCreate(ls[sbi]->get_solution(), is, output_vec[sbi], PETSC_NULL, &output_scatter);
382  VecScatterBegin(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
383  VecScatterEnd(output_scatter, ls[sbi]->get_solution(), output_vec[sbi], INSERT_VALUES, SCATTER_FORWARD);
384  VecScatterDestroy(&(output_scatter));
385  ISDestroy(&(is));
386  }
387 }
388 
389 
390 
391 template<class Model>
393 {
394  data_.set_time(*time_);
395 
396  // set initial conditions
397  set_initial_condition();
398  for (unsigned int sbi = 0; sbi < n_subst_; sbi++)
399  ( (LinSys_PETSC *)ls[sbi] )->set_initial_guess_nonzero();
400 
401  output_data();
402 }
403 
404 
405 template<class Model>
407 {
408  START_TIMER("DG-ONE STEP");
409 
410  time_->next_time();
411  time_->view("TDG");
412 
413  START_TIMER("data reinit");
414  data_.set_time(*time_);
415  END_TIMER("data reinit");
416 
417  // check first time assembly - needs preallocation
418  if (!allocation_done)
419  {
420  // preallocate mass matrix
421  ls_dt->start_allocation();
422  assemble_mass_matrix();
423  mass_matrix = NULL;
424 
425  // preallocate system matrix
426  for (unsigned int i=0; i<n_subst_; i++)
427  {
428  ls[i]->start_allocation();
429  stiffness_matrix[i] = NULL;
430  rhs[i] = NULL;
431  }
432  assemble_stiffness_matrix();
433  set_sources();
434  set_boundary_conditions();
435 
436  allocation_done = true;
437  }
438 
439  // assemble mass matrix
440  if (mass_matrix == NULL || data_.subset(FieldFlag::in_time_term).changed() )
441  {
442  ls_dt->start_add_assembly();
443  ls_dt->mat_zero_entries();
444  assemble_mass_matrix();
445  ls_dt->finish_assembly();
446  mass_matrix = *(ls_dt->get_matrix());
447  }
448 
449  // assemble stiffness matrix
450  if (stiffness_matrix[0] == NULL
451  || data_.subset(FieldFlag::in_main_matrix).changed()
452  || Model::flux_changed)
453  {
454  // new fluxes can change the location of Neumann boundary,
455  // thus stiffness matrix must be reassembled
456  for (unsigned int i=0; i<n_subst_; i++)
457  {
458  ls[i]->start_add_assembly();
459  ls[i]->mat_zero_entries();
460  }
461  assemble_stiffness_matrix();
462  for (unsigned int i=0; i<n_subst_; i++)
463  {
464  ls[i]->finish_assembly();
465 
466  if (stiffness_matrix[i] == NULL)
467  MatConvert(*( ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
468  else
469  MatCopy(*( ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
470  }
471  }
472 
473  // assemble right hand side (due to sources and boundary conditions)
474  if (rhs[0] == NULL
475  || data_.subset(FieldFlag::in_rhs).changed()
476  || Model::flux_changed)
477  {
478  for (unsigned int i=0; i<n_subst_; i++)
479  {
480  ls[i]->start_add_assembly();
481  ls[i]->rhs_zero_entries();
482  }
483  set_sources();
484  set_boundary_conditions();
485  for (unsigned int i=0; i<n_subst_; i++)
486  {
487  ls[i]->finish_assembly();
488 
489  VecDuplicate(*( ls[i]->get_rhs() ), &rhs[i]);
490  VecCopy(*( ls[i]->get_rhs() ), rhs[i]);
491  }
492  }
493 
494  Model::flux_changed = false;
495 
496 
497  /* Apply backward Euler time integration.
498  *
499  * Denoting A the stiffness matrix and M the mass matrix, the algebraic system at the k-th time level reads
500  *
501  * (1/dt M + A)u^k = f + 1/dt M.u^{k-1}
502  *
503  * Hence we modify at each time level the right hand side:
504  *
505  * f^k = f + 1/dt M u^{k-1},
506  *
507  * where f stands for the term stemming from the force and boundary conditions.
508  * Accordingly, we set
509  *
510  * A^k = A + 1/dt M.
511  *
512  */
513  Mat m;
514  START_TIMER("solve");
515  for (unsigned int i=0; i<n_subst_; i++)
516  {
517  MatConvert(stiffness_matrix[i], MATSAME, MAT_INITIAL_MATRIX, &m);
518  MatAXPY(m, 1./time_->dt(), mass_matrix, SUBSET_NONZERO_PATTERN);
519  ls[i]->set_matrix(m, DIFFERENT_NONZERO_PATTERN);
520  Vec y,w;
521  VecDuplicate(rhs[i], &y);
522  VecDuplicate(rhs[i], &w);
523  MatMult(mass_matrix, ls[i]->get_solution(), y);
524  VecWAXPY(w, 1./time_->dt(), y, rhs[i]);
525  ls[i]->set_rhs(w);
526 
527  VecDestroy(&y);
528  VecDestroy(&w);
529  MatDestroy(&m);
530 
531  ls[i]->solve();
532  }
533  END_TIMER("solve");
534 
535  if (mass_balance() != NULL)
536  mass_balance()->calculate(time_->t());
537 
538  END_TIMER("DG-ONE STEP");
539 }
540 
541 
542 template<class Model>
544 {
545  // So far the velocity_vector contains zeros, so we ignore it.
546  // Instead we use the value Side.flux.
547 
548  mh_dh = &dh;
549  Model::flux_changed = true;
550 
551 }
552 
553 
554 template<class Model>
556 {
557  if (!time_->is_current( time_->marks().type_output() )) return;
558 
559  START_TIMER("DG-OUTPUT");
560 
561  // gather the solution from all processors
562  output_vector_gather();
563  data_.subset(FieldFlag::allow_output).set_time( *time_);
564  data_.output(output_stream);
565  output_stream->write_time_frame();
566 
567  if (mass_balance() != NULL)
568  mass_balance()->output(time_->t());
569 
570  END_TIMER("DG-OUTPUT");
571 }
572 
573 
574 template<class Model>
576 {
577  START_TIMER("assemble_mass");
578  assemble_mass_matrix<1>();
579  assemble_mass_matrix<2>();
580  assemble_mass_matrix<3>();
581  END_TIMER("assemble_mass");
582 }
583 
584 
585 template<class Model> template<unsigned int dim>
587 {
588  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(), update_values | update_JxW_values | update_quadrature_points);
589  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
590  unsigned int dof_indices[ndofs];
591  PetscScalar local_mass_matrix[ndofs*ndofs];
592 
593  // assemble integral over elements
594  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
595  {
596  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
597  if (cell->dim() != dim) continue;
598 
599  fe_values.reinit(cell);
600  feo->dh()->get_dof_indices(cell, dof_indices);
601  ElementAccessor<3> ele_acc = cell->element_accessor();
602 
603  Model::compute_mass_matrix_coefficient(fe_values.point_list(), ele_acc, mm_coef);
604 
605  // assemble the local mass matrix
606  for (unsigned int i=0; i<ndofs; i++)
607  {
608  for (unsigned int j=0; j<ndofs; j++)
609  {
610  local_mass_matrix[i*ndofs+j] = 0;
611  for (unsigned int k=0; k<qsize; k++)
612  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);
613  }
614  }
615 
616  ls_dt->mat_set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, local_mass_matrix);
617  }
618 }
619 
620 
621 
622 
623 template<class Model>
625 {
626  START_TIMER("assemble_stiffness");
627  START_TIMER("assemble_volume_integrals");
628  assemble_volume_integrals<1>();
629  assemble_volume_integrals<2>();
630  assemble_volume_integrals<3>();
631  END_TIMER("assemble_volume_integrals");
632 
633  START_TIMER("assemble_fluxes_boundary");
634  assemble_fluxes_boundary<1>();
635  assemble_fluxes_boundary<2>();
636  assemble_fluxes_boundary<3>();
637  END_TIMER("assemble_fluxes_boundary");
638 
639  START_TIMER("assemble_fluxes_elem_elem");
640  assemble_fluxes_element_element<1>();
641  assemble_fluxes_element_element<2>();
642  assemble_fluxes_element_element<3>();
643  END_TIMER("assemble_fluxes_elem_elem");
644 
645  START_TIMER("assemble_fluxes_elem_side");
646  assemble_fluxes_element_side<1>();
647  assemble_fluxes_element_side<2>();
648  assemble_fluxes_element_side<3>();
649  END_TIMER("assemble_fluxes_elem_side");
650  END_TIMER("assemble_stiffness");
651 }
652 
653 
654 
655 template<class Model>
656 template<unsigned int dim>
658 {
659  FEValues<dim,3> fv_rt(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe_rt<dim>(),
661  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
663  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
664  unsigned int dof_indices[ndofs];
665  vector<arma::vec3> velocity(qsize);
666  vector<arma::vec> sources_sigma(qsize);
667  PetscScalar local_matrix[ndofs*ndofs];
668 
669  // assemble integral over elements
670  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
671  {
672  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
673  if (cell->dim() != dim) continue;
674 
675  fe_values.reinit(cell);
676  fv_rt.reinit(cell);
677  ElementAccessor<3> ele_acc = cell->element_accessor();
678  feo->dh()->get_dof_indices(cell, dof_indices);
679 
680  calculate_velocity(cell, velocity, fv_rt);
681  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), velocity, ele_acc, ad_coef, dif_coef);
682  Model::compute_sources_sigma(fe_values.point_list(), ele_acc, sources_sigma);
683 
684  // assemble the local stiffness matrix
685  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
686  {
687  for (unsigned int i=0; i<ndofs; i++)
688  for (unsigned int j=0; j<ndofs; j++)
689  local_matrix[i*ndofs+j] = 0;
690 
691  for (unsigned int k=0; k<qsize; k++)
692  {
693  for (unsigned int i=0; i<ndofs; i++)
694  {
695  arma::vec3 Kt_grad_i = dif_coef[sbi][k].t()*fe_values.shape_grad(i,k);
696  double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.shape_grad(i,k));
697 
698  for (unsigned int j=0; j<ndofs; j++)
699  local_matrix[i*ndofs+j] += (arma::dot(Kt_grad_i, fe_values.shape_grad(j,k))
700  -fe_values.shape_value(j,k)*ad_dot_grad_i
701  +sources_sigma[k][sbi]*fe_values.shape_value(j,k)*fe_values.shape_value(i,k))*fe_values.JxW(k);
702  }
703  }
704  ls[sbi]->mat_set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, local_matrix);
705  }
706  }
707 }
708 
709 
710 template<class Model>
712 {
713  START_TIMER("assemble_sources");
714  set_sources<1>();
715  set_sources<2>();
716  set_sources<3>();
717  END_TIMER("assemble_sources");
718 }
719 
720 template<class Model>
721 template<unsigned int dim>
723 {
724  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
726  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
727  vector<arma::vec> sources_conc(qsize), sources_density(qsize), sources_sigma(qsize);
728  unsigned int dof_indices[ndofs];
729  PetscScalar local_rhs[ndofs];
730  double source;
731 
732  // assemble integral over elements
733  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
734  {
735  typename DOFHandlerBase::CellIterator cell = mesh_->element(feo->dh()->el_index(i_cell));
736  if (cell->dim() != dim) continue;
737 
738  fe_values.reinit(cell);
739  feo->dh()->get_dof_indices(cell, dof_indices);
740 
741  Model::compute_source_coefficients(fe_values.point_list(), cell->element_accessor(), sources_conc, sources_density, sources_sigma);
742 
743  // assemble the local stiffness matrix
744  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
745  {
746  for (unsigned int i=0; i<ndofs; i++)
747  local_rhs[i] = 0;
748 
749  // compute sources
750  for (unsigned int k=0; k<qsize; k++)
751  {
752  source = (sources_density[k][sbi] + sources_conc[k][sbi]*sources_sigma[k][sbi])*fe_values.JxW(k);
753 
754  for (unsigned int i=0; i<ndofs; i++)
755  local_rhs[i] += source*fe_values.shape_value(i,k);
756  }
757  ls[sbi]->rhs_set_values(ndofs, (int *)dof_indices, local_rhs);
758  }
759  }
760 }
761 
762 
763 
764 template<class Model>
765 template<unsigned int dim>
767 {
768  vector<FESideValues<dim,3>*> fe_values;
769  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
770  update_values);
771  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size(),
772  n_max_sides = ad_coef_edg.size();
773  vector<unsigned int*> side_dof_indices;
774  PetscScalar local_matrix[ndofs*ndofs];
775  vector<vector<arma::vec3> > side_velocity(n_max_sides);
776  vector<arma::vec> dg_penalty(n_max_sides);
777  double gamma_l, omega[2], transport_flux;
778 
779  for (unsigned int sid=0; sid<n_max_sides; sid++)
780  {
781  side_dof_indices.push_back(new unsigned int[ndofs]);
782  fe_values.push_back(new FESideValues<dim,3>(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
784  }
785 
786  // assemble integral over sides
787  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
788  {
789  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
790  if (edg->n_sides < 2 || edg->side(0)->element()->dim() != dim) continue;
791 
792  for (int sid=0; sid<edg->n_sides; sid++)
793  {
794  typename DOFHandlerBase::CellIterator cell = edg->side(sid)->element();
795  ElementAccessor<3> ele_acc = cell->element_accessor();
796  feo->dh()->get_dof_indices(cell, side_dof_indices[sid]);
797  fe_values[sid]->reinit(cell, edg->side(sid)->el_idx());
798  fsv_rt.reinit(cell, edg->side(sid)->el_idx());
799  calculate_velocity(cell, side_velocity[sid], fsv_rt);
800  Model::compute_advection_diffusion_coefficients(fe_values[sid]->point_list(), side_velocity[sid], ele_acc, ad_coef_edg[sid], dif_coef_edg[sid]);
801  dg_penalty[sid] = data_.dg_penalty.value(cell->centre(), ele_acc);
802  }
803 
804  // fluxes and penalty
805  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
806  {
807  vector<double> fluxes(edg->n_sides);
808  for (int sid=0; sid<edg->n_sides; sid++)
809  {
810  fluxes[sid] = 0;
811  for (unsigned int k=0; k<qsize; k++)
812  fluxes[sid] += arma::dot(ad_coef_edg[sid][sbi][k], fe_values[sid]->normal_vector(k))*fe_values[sid]->JxW(k);
813  fluxes[sid] /= edg->side(sid)->measure();
814  }
815 
816  for (int s1=0; s1<edg->n_sides; s1++)
817  {
818  for (int s2=s1+1; s2<edg->n_sides; s2++)
819  {
820  ASSERT(edg->side(s1)->valid(), "Invalid side of edge.");
821  if (!feo->dh()->el_is_local(edg->side(s1)->element().index())
822  && !feo->dh()->el_is_local(edg->side(s2)->element().index())) continue;
823 
824  arma::vec3 nv = fe_values[s1]->normal_vector(0);
825 
826  // set up the parameters for DG method
827  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);
828 
829  int sd[2];
830  sd[0] = s1;
831  sd[1] = s2;
832 
833 #define AVERAGE(i,k,side_id) (fe_values[sd[side_id]]->shape_value(i,k)*0.5)
834 #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])
835 #define JUMP(i,k,side_id) ((side_id==0?1:-1)*fe_values[sd[side_id]]->shape_value(i,k))
836 
837  // For selected pair of elements:
838  for (int n=0; n<2; n++)
839  {
840  if (!feo->dh()->el_is_local(edg->side(sd[n])->element().index())) continue;
841 
842  for (int m=0; m<2; m++)
843  {
844  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
845  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
846  local_matrix[i*fe_values[sd[m]]->n_dofs()+j] = 0;
847 
848  for (unsigned int k=0; k<qsize; k++)
849  {
850  double flux_times_JxW = transport_flux*fe_values[0]->JxW(k);
851  double gamma_times_JxW = gamma_l*fe_values[0]->JxW(k);
852 
853  for (unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
854  {
855  double flux_JxW_jump_i = flux_times_JxW*JUMP(i,k,n);
856  double gamma_JxW_jump_i = gamma_times_JxW*JUMP(i,k,n);
857  double JxW_jump_i = fe_values[0]->JxW(k)*JUMP(i,k,n);
858  double JxW_var_wavg_i = fe_values[0]->JxW(k)*WAVERAGE(i,k,n)*dg_variant;
859 
860  for (unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
861  {
862  int index = i*fe_values[sd[m]]->n_dofs()+j;
863 
864  // flux due to transport (applied on interior edges) (average times jump)
865  local_matrix[index] += flux_JxW_jump_i*AVERAGE(j,k,m);
866 
867  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
868  local_matrix[index] += gamma_JxW_jump_i*JUMP(j,k,m);
869 
870  // terms due to diffusion
871  local_matrix[index] -= WAVERAGE(j,k,m)*JxW_jump_i;
872  local_matrix[index] -= JUMP(j,k,m)*JxW_var_wavg_i;
873  }
874  }
875  }
876  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);
877  }
878  }
879 #undef AVERAGE
880 #undef WAVERAGE
881 #undef JUMP
882  }
883  }
884  }
885  }
886 
887  for (unsigned int i=0; i<n_max_sides; i++)
888  {
889  delete fe_values[i];
890  delete[] side_dof_indices[i];
891  }
892 }
893 
894 
895 template<class Model>
896 template<unsigned int dim>
898 {
899  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
901  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
902  update_values);
903  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
904  unsigned int side_dof_indices[ndofs];
905  PetscScalar local_matrix[ndofs*ndofs];
906  vector<arma::vec3> side_velocity;
907  vector<arma::vec> robin_sigma(qsize);
908  arma::vec dg_penalty;
909  double gamma_l;
910 
911  // assemble boundary integral
912  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
913  {
914  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
915  if (edg->n_sides > 1) continue;
916  // check spatial dimension
917  if (edg->side(0)->dim() != dim-1) continue;
918  // skip edges lying not on the boundary
919  if (edg->side(0)->cond() == NULL) continue;
920 
921  SideIter side = edg->side(0);
922  typename DOFHandlerBase::CellIterator cell = side->element();
923  ElementAccessor<3> ele_acc = cell->element_accessor();
924  feo->dh()->get_dof_indices(cell, side_dof_indices);
925  fe_values_side.reinit(cell, side->el_idx());
926  fsv_rt.reinit(cell, side->el_idx());
927 
928  calculate_velocity(cell, side_velocity, fsv_rt);
929  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), side_velocity, ele_acc, ad_coef, dif_coef);
930  dg_penalty = data_.dg_penalty.value(cell->centre(), ele_acc);
931  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), side->cond()->element_accessor());
932  data_.bc_robin_sigma.value_list(fe_values_side.point_list(), side->cond()->element_accessor(), robin_sigma);
933 
934  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
935  {
936  for (unsigned int i=0; i<ndofs; i++)
937  for (unsigned int j=0; j<ndofs; j++)
938  local_matrix[i*ndofs+j] = 0;
939 
940  // On Neumann boundaries we have only term from integrating by parts the advective term,
941  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
942  double side_flux = 0;
943  for (unsigned int k=0; k<qsize; k++)
944  side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.normal_vector(k))*fe_values_side.JxW(k);
945  double transport_flux = side_flux/side->measure();
946 
947  if ((bc_type[sbi] == EqData::dirichlet) || (bc_type[sbi] == EqData::inflow ))
948  {
949  // set up the parameters for DG method
950  set_DG_parameters_boundary(side, qsize, dif_coef[sbi], transport_flux, fe_values_side.normal_vector(0), dg_penalty[sbi], gamma_l);
951  gamma[sbi][side->cond_idx()] = gamma_l;
952  if (bc_type[sbi] == EqData::dirichlet || side_flux < -mh_dh->precision())
953  transport_flux += gamma_l;
954  }
955 
956  // fluxes and penalty
957  for (unsigned int k=0; k<qsize; k++)
958  {
959  double flux_times_JxW;
960  if (bc_type[sbi] == EqData::robin)
961  flux_times_JxW = (transport_flux + robin_sigma[k][sbi])*fe_values_side.JxW(k);
962  else
963  flux_times_JxW = transport_flux*fe_values_side.JxW(k);
964 
965  for (unsigned int i=0; i<ndofs; i++)
966  for (unsigned int j=0; j<ndofs; j++)
967  {
968  // flux due to advection and penalty
969  local_matrix[i*ndofs+j] += flux_times_JxW*fe_values_side.shape_value(i,k)*fe_values_side.shape_value(j,k);
970 
971  // flux due to diffusion (only on dirichlet and inflow boundary)
972  if (bc_type[sbi] == EqData::dirichlet || (bc_type[sbi] == EqData::inflow && side_flux < -mh_dh->precision()))
973  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)
974  + 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
975  )*fe_values_side.JxW(k);
976  }
977  }
978  ls[sbi]->mat_set_values(ndofs, (int *)side_dof_indices, ndofs, (int *)side_dof_indices, local_matrix);
979  }
980  }
981 }
982 
983 
984 template<class Model>
985 template<unsigned int dim>
987 {
988 
989  if (dim == 1) return;
990  FEValues<dim-1,3> fe_values_vb(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe<dim-1>(),
992  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
994  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
995  update_values);
996  FEValues<dim-1,3> fv_rt(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe_rt<dim-1>(),
997  update_values);
998 
999  vector<FEValuesSpaceBase<3>*> fv_sb(2);
1000  const unsigned int ndofs = feo->fe<dim>()->n_dofs(); // number of local dofs
1001  const unsigned int qsize = feo->q<dim-1>()->size(); // number of quadrature points
1002  unsigned int side_dof_indices[2*ndofs], n_dofs[2];
1003  vector<arma::vec3> velocity_higher, velocity_lower;
1004  vector<arma::vec> frac_sigma(qsize);
1005  vector<double> csection_lower(qsize), csection_higher(qsize), mm_coef_lower(qsize), mm_coef_higher(qsize);
1006  PetscScalar local_matrix[4*ndofs*ndofs];
1007  double comm_flux[2][2];
1008 
1009  // index 0 = element with lower dimension,
1010  // index 1 = side of element with higher dimension
1011  fv_sb[0] = &fe_values_vb;
1012  fv_sb[1] = &fe_values_side;
1013 
1014  // assemble integral over sides
1015  for (unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
1016  {
1017  Neighbour *nb = &mesh_->vb_neighbours_[feo->dh()->nb_index(inb)];
1018  // skip neighbours of different dimension
1019  if (nb->element()->dim() != dim-1) continue;
1020 
1021  typename DOFHandlerBase::CellIterator cell_sub = mesh_->element.full_iter(nb->element());
1022  feo->dh()->get_dof_indices(cell_sub, side_dof_indices);
1023  fe_values_vb.reinit(cell_sub);
1024  n_dofs[0] = fv_sb[0]->n_dofs();
1025 
1026  typename DOFHandlerBase::CellIterator cell = nb->side()->element();
1027  feo->dh()->get_dof_indices(cell, side_dof_indices+n_dofs[0]);
1028  fe_values_side.reinit(cell, nb->side()->el_idx());
1029  n_dofs[1] = fv_sb[1]->n_dofs();
1030 
1031  // Element id's for testing if they belong to local partition.
1032  int element_id[2];
1033  element_id[0] = cell_sub.index();
1034  element_id[1] = cell.index();
1035 
1036  fsv_rt.reinit(cell, nb->side()->el_idx());
1037  fv_rt.reinit(cell_sub);
1038  calculate_velocity(cell, velocity_higher, fsv_rt);
1039  calculate_velocity(cell_sub, velocity_lower, fv_rt);
1040  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_lower, cell_sub->element_accessor(), ad_coef_edg[0], dif_coef_edg[0]);
1041  Model::compute_advection_diffusion_coefficients(fe_values_vb.point_list(), velocity_higher, cell->element_accessor(), ad_coef_edg[1], dif_coef_edg[1]);
1042  Model::compute_mass_matrix_coefficient(fe_values_vb.point_list(), cell_sub->element_accessor(), mm_coef_lower);
1043  Model::compute_mass_matrix_coefficient(fe_values_vb.point_list(), cell->element_accessor(), mm_coef_higher);
1044  data_.cross_section.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), csection_lower);
1045  data_.cross_section.value_list(fe_values_vb.point_list(), cell->element_accessor(), csection_higher);
1046  data_.fracture_sigma.value_list(fe_values_vb.point_list(), cell_sub->element_accessor(), frac_sigma);
1047 
1048  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1049  {
1050  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
1051  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
1052  local_matrix[i*(n_dofs[0]+n_dofs[1])+j] = 0;
1053 
1054  // set transmission conditions
1055  for (unsigned int k=0; k<qsize; k++)
1056  {
1057  /* The communication flux has two parts:
1058  * - "diffusive" term containing sigma
1059  * - "advective" term representing usual upwind
1060  *
1061  * The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
1062  * than b and A in the manual.
1063  * In calculation of sigma there appears one more csection_lower in the denominator.
1064  */
1065  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))*
1066  2*csection_higher[k]*csection_higher[k]/(csection_lower[k]*csection_lower[k]);
1067 
1068  // Since mm_coef_* contains cross section, we have to divide by it.
1069  double transport_flux = arma::dot(ad_coef_edg[1][sbi][k], fe_values_side.normal_vector(k));
1070  double por_lower_over_higher = mm_coef_lower[k]*csection_higher[k]/(mm_coef_higher[k]*csection_lower[k]);
1071 
1072  comm_flux[0][0] = (sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1073  comm_flux[0][1] = -(sigma-min(0.,transport_flux*por_lower_over_higher))*fv_sb[0]->JxW(k);
1074  comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1075  comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1076 
1077  for (int n=0; n<2; n++)
1078  {
1079  if (!feo->dh()->el_is_local(element_id[n])) continue;
1080 
1081  for (unsigned int i=0; i<n_dofs[n]; i++)
1082  for (int m=0; m<2; m++)
1083  for (unsigned int j=0; j<n_dofs[m]; j++)
1084  local_matrix[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
1085  comm_flux[m][n]*fv_sb[m]->shape_value(j,k)*fv_sb[n]->shape_value(i,k);
1086  }
1087  }
1088  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);
1089  }
1090  }
1091 
1092 }
1093 
1094 
1095 
1096 
1097 
1098 
1099 template<class Model>
1101 {
1102  START_TIMER("assemble_bc");
1103  set_boundary_conditions<1>();
1104  set_boundary_conditions<2>();
1105  set_boundary_conditions<3>();
1106  END_TIMER("assemble_bc");
1107 }
1108 
1109 
1110 template<class Model>
1111 template<unsigned int dim>
1113 {
1114  FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1116  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1117  update_values);
1118  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1119  unsigned int side_dof_indices[ndofs];
1120  double local_rhs[ndofs];
1121  vector<arma::vec> bc_values(qsize), bc_fluxes(qsize), bc_sigma(qsize);
1122  vector<arma::vec3> velocity;
1123 
1124  for (unsigned int i=0; i<qsize; i++)
1125  {
1126  bc_values[i].resize(n_subst_);
1127  bc_fluxes[i].resize(n_subst_);
1128  bc_sigma[i].resize(n_subst_);
1129  }
1130 
1131  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1132  {
1133  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
1134  if (edg->n_sides > 1) continue;
1135  if (edg->side(0)->dim() != dim-1) continue;
1136  // skip edges lying not on the boundary
1137  if (edg->side(0)->cond() == NULL) continue;
1138 
1139  SideIter side = edg->side(0);
1140  typename DOFHandlerBase::CellIterator cell = mesh().element.full_iter(side->element());
1141  ElementAccessor<3> ele_acc = side->cond()->element_accessor();
1142 
1143  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), ele_acc);
1144 
1145  fe_values_side.reinit(cell, side->el_idx());
1146  fsv_rt.reinit(cell, side->el_idx());
1147  calculate_velocity(cell, velocity, fsv_rt);
1148 
1149  Model::compute_advection_diffusion_coefficients(fe_values_side.point_list(), velocity, side->element()->element_accessor(), ad_coef, dif_coef);
1150  Model::compute_dirichlet_bc(fe_values_side.point_list(), ele_acc, bc_values);
1151  data_.bc_flux.value_list(fe_values_side.point_list(), ele_acc, bc_fluxes);
1152  data_.bc_robin_sigma.value_list(fe_values_side.point_list(), ele_acc, bc_sigma);
1153 
1154  feo->dh()->get_dof_indices(cell, side_dof_indices);
1155 
1156  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1157  {
1158  for (unsigned int i=0; i<ndofs; i++) local_rhs[i] = 0;
1159 
1160  for (unsigned int k=0; k<qsize; k++)
1161  {
1162  double bc_term = 0;
1163  arma::vec3 bc_grad;
1164  bc_grad.zeros();
1165  if ((bc_type[sbi] == EqData::inflow && mh_dh->side_flux( *edg->side(0) ) < -mh_dh->precision())
1166  || (bc_type[sbi] == EqData::dirichlet))
1167  {
1168  bc_term = gamma[sbi][side->cond_idx()]*bc_values[k][sbi]*fe_values_side.JxW(k);
1169  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));
1170  }
1171  else if (bc_type[sbi] == EqData::neumann)
1172  {
1173  bc_term = bc_fluxes[k][sbi]*fe_values_side.JxW(k);
1174  }
1175  else if (bc_type[sbi] == EqData::robin)
1176  {
1177  bc_term = bc_sigma[k][sbi]*bc_values[k][sbi]*fe_values_side.JxW(k);
1178  }
1179 
1180  for (unsigned int i=0; i<ndofs; i++)
1181  local_rhs[i] += bc_term*fe_values_side.shape_value(i,k)
1182  + arma::dot(bc_grad,fe_values_side.shape_grad(i,k));
1183  }
1184  ls[sbi]->rhs_set_values(ndofs, (int *)side_dof_indices, local_rhs);
1185  }
1186  }
1187 }
1188 
1189 
1190 
1191 // TODO: The detection of side number from SideIter
1192 // in TransportDG::calculate_velocity()
1193 // should be done with help of class RefElement. This however requires
1194 // that the MH_DofHandler uses the node/side ordering defined in
1195 // the respective RefElement.
1196 template<class Model>
1197 template<unsigned int dim>
1199 {
1200  std::map<const Node*, int> node_nums;
1201  for (unsigned int i=0; i<cell->n_nodes(); i++)
1202  node_nums[cell->node[i]] = i;
1203 
1204  velocity.resize(fv.n_points());
1205 
1206  for (unsigned int k=0; k<fv.n_points(); k++)
1207  {
1208  velocity[k].zeros();
1209  for (unsigned int sid=0; sid<cell->n_sides(); sid++)
1210  {
1211  if (cell->side(sid)->dim() != dim-1) continue;
1212  int num = dim*(dim+1)/2;
1213  for (unsigned int i=0; i<cell->side(sid)->n_nodes(); i++)
1214  num -= node_nums[cell->side(sid)->node(i)];
1215  velocity[k] += fv.shape_vector(num,k) * mh_dh->side_flux( *(cell->side(sid)) );
1216  }
1217  }
1218 }
1219 
1220 
1221 
1222 
1223 
1224 template<class Model>
1226  const int s1,
1227  const int s2,
1228  const int K_size,
1229  const vector<arma::mat33> &K1,
1230  const vector<arma::mat33> &K2,
1231  const vector<double> &fluxes,
1232  const arma::vec3 &normal_vector,
1233  const double alpha1,
1234  const double alpha2,
1235  double &gamma,
1236  double *omega,
1237  double &transport_flux)
1238 {
1239  double delta[2];
1240  double h = 0;
1241  double local_alpha = 0;
1242 
1243  ASSERT(edg.side(s1)->valid(), "Invalid side of an edge.");
1244  SideIter s = edg.side(s1);
1245 
1246  // calculate the side diameter
1247  if (s->dim() == 0)
1248  {
1249  h = 1;
1250  }
1251  else
1252  {
1253  for (unsigned int i=0; i<s->n_nodes(); i++)
1254  for (unsigned int j=i+1; j<s->n_nodes(); j++)
1255  h = max(h, s->node(i)->distance(*s->node(j)));
1256  }
1257 
1258  // calculate the total in- and out-flux through the edge
1259  double pflux = 0, nflux = 0;
1260  for (int i=0; i<edg.n_sides; i++)
1261  {
1262  if (fluxes[i] > 0)
1263  pflux += fluxes[i];
1264  else
1265  nflux += fluxes[i];
1266  }
1267 
1268  // calculate the flux from s1 to s2
1269  if (fluxes[s2] > 0 && fluxes[s1] < 0 && s1 < s2)
1270  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
1271  else if (fluxes[s2] < 0 && fluxes[s1] > 0 && s1 < s2)
1272  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
1273  else if (s1==s2)
1274  transport_flux = fluxes[s1];
1275  else
1276  transport_flux = 0;
1277 
1278  gamma = 0.5*fabs(transport_flux);
1279 
1280 
1281  // determine local DG penalty
1282  local_alpha = max(alpha1, alpha2);
1283 
1284  if (s1 == s2)
1285  {
1286  omega[0] = 1;
1287 
1288  // delta is set to the average value of Kn.n on the side
1289  delta[0] = 0;
1290  for (int k=0; k<K_size; k++)
1291  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1292  delta[0] /= K_size;
1293 
1294  gamma += local_alpha/h*(delta[0]);
1295  }
1296  else
1297  {
1298  delta[0] = 0;
1299  delta[1] = 0;
1300  for (int k=0; k<K_size; k++)
1301  {
1302  delta[0] += dot(K1[k]*normal_vector,normal_vector);
1303  delta[1] += dot(K2[k]*normal_vector,normal_vector);
1304  }
1305  delta[0] /= K_size;
1306  delta[1] /= K_size;
1307 
1308  double delta_sum = delta[0] + delta[1];
1309 
1310  if (delta_sum > numeric_limits<double>::epsilon())
1311  {
1312  omega[0] = delta[1]/delta_sum;
1313  omega[1] = delta[0]/delta_sum;
1314  gamma += local_alpha/h*(delta[0]*delta[1]/delta_sum);
1315  }
1316  else
1317  for (int i=0; i<2; i++) omega[i] = 0;
1318  }
1319 }
1320 
1321 
1322 
1323 
1324 
1325 
1326 template<class Model>
1328  const int K_size,
1329  const vector<arma::mat33> &K,
1330  const double flux,
1331  const arma::vec3 &normal_vector,
1332  const double alpha,
1333  double &gamma)
1334 {
1335  double delta = 0, h = 0;
1336 
1337  // calculate the side diameter
1338  if (side->dim() == 0)
1339  {
1340  h = 1;
1341  }
1342  else
1343  {
1344  for (unsigned int i=0; i<side->n_nodes(); i++)
1345  for (unsigned int j=i+1; j<side->n_nodes(); j++)
1346  h = max(h, side->node(i)->distance( *side->node(j) ));
1347  }
1348 
1349  // delta is set to the average value of Kn.n on the side
1350  for (int k=0; k<K_size; k++)
1351  delta += dot(K[k]*normal_vector,normal_vector);
1352  delta /= K_size;
1353 
1354  gamma = 0.5*fabs(flux) + alpha/h*delta;
1355 }
1356 
1357 
1358 
1359 
1360 
1361 template<class Model>
1363 {
1364  START_TIMER("set_init_cond");
1365  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1366  ls[sbi]->start_allocation();
1367  prepare_initial_condition<1>();
1368  prepare_initial_condition<2>();
1369  prepare_initial_condition<3>();
1370 
1371  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1372  ls[sbi]->start_add_assembly();
1373  prepare_initial_condition<1>();
1374  prepare_initial_condition<2>();
1375  prepare_initial_condition<3>();
1376 
1377  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1378  {
1379  ls[sbi]->finish_assembly();
1380  ls[sbi]->solve();
1381  }
1382  END_TIMER("set_init_cond");
1383 }
1384 
1385 template<class Model>
1386 template<unsigned int dim>
1388 {
1389  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1391  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1392  unsigned int dof_indices[ndofs];
1393  double matrix[ndofs*ndofs], rhs[ndofs];
1394  std::vector<arma::vec> init_values(qsize);
1395 
1396  for (unsigned int k=0; k<qsize; k++)
1397  init_values[k].resize(n_subst_);
1398 
1399  for (unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
1400  {
1401  typename DOFHandlerBase::CellIterator elem = mesh_->element(feo->dh()->el_index(i_cell));
1402  if (elem->dim() != dim) continue;
1403 
1404  ElementAccessor<3> ele_acc = elem->element_accessor();
1405  feo->dh()->get_dof_indices(elem, dof_indices);
1406  fe_values.reinit(elem);
1407 
1408  Model::compute_init_cond(fe_values.point_list(), ele_acc, init_values);
1409 
1410  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1411  {
1412  for (unsigned int i=0; i<ndofs; i++)
1413  {
1414  rhs[i] = 0;
1415  for (unsigned int j=0; j<ndofs; j++)
1416  matrix[i*ndofs+j] = 0;
1417  }
1418 
1419  for (unsigned int k=0; k<qsize; k++)
1420  {
1421  double rhs_term = init_values[k](sbi)*fe_values.JxW(k);
1422 
1423  for (unsigned int i=0; i<ndofs; i++)
1424  {
1425  for (unsigned int j=0; j<ndofs; j++)
1426  matrix[i*ndofs+j] += fe_values.shape_value(i,k)*fe_values.shape_value(j,k)*fe_values.JxW(k);
1427 
1428  rhs[i] += fe_values.shape_value(i,k)*rhs_term;
1429  }
1430  }
1431  ls[sbi]->set_values(ndofs, (int *)dof_indices, ndofs, (int *)dof_indices, matrix, rhs);
1432  }
1433  }
1434 }
1435 
1436 
1437 template<class Model>
1438 void TransportDG<Model>::calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance)
1439 {
1440  calc_fluxes<1>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1441  calc_fluxes<2>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1442  calc_fluxes<3>(bcd_balance, bcd_plus_balance, bcd_minus_balance);
1443 }
1444 
1445 template<class Model>
1446 template<unsigned int dim>
1447 void TransportDG<Model>::calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance)
1448 {
1449  FESideValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
1451  FESideValues<dim,3> fsv_rt(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe_rt<dim>(),
1452  update_values);
1453  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1454  unsigned int dof_indices[ndofs];
1455  vector<arma::vec3> side_velocity(qsize);
1456  vector<arma::vec> bc_values(qsize);
1457  arma::vec3 conc_grad;
1458 
1459  for (unsigned int k=0; k<qsize; k++)
1460  bc_values[k].resize(n_subst_);
1461 
1462  for (unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1463  {
1464  Edge *edg = &mesh_->edges[feo->dh()->edge_index(iedg)];
1465  if (edg->n_sides > 1) continue;
1466  if (edg->side(0)->dim() != dim-1) continue;
1467  // skip edges lying not on the boundary
1468  if (edg->side(0)->cond() == NULL) continue;
1469 
1470  SideIter side = edg->side(0);
1472  ElementAccessor<3> ele_acc = cell->element_accessor();
1473  Region r = side->cond()->element_accessor().region();
1474  if (! r.is_valid()) xprintf(Msg, "Invalid region, ele % d, edg: % d\n", side->cond()->bc_ele_idx_, side->cond()->edge_idx_);
1475  unsigned int bc_region_idx = r.boundary_idx();
1476 
1477  fe_values.reinit(cell, side->el_idx());
1478  fsv_rt.reinit(cell, side->el_idx());
1479  feo->dh()->get_dof_indices(cell, dof_indices);
1480  calculate_velocity(cell, side_velocity, fsv_rt);
1481  Model::compute_advection_diffusion_coefficients(fe_values.point_list(), side_velocity, ele_acc, ad_coef, dif_coef);
1482  Model::compute_dirichlet_bc(fe_values.point_list(), side->cond()->element_accessor(), bc_values);
1483  arma::uvec bc_type = data_.bc_type.value(side->cond()->element()->centre(), side->cond()->element_accessor());
1484 
1485  for (unsigned int sbi=0; sbi<n_subst_; sbi++)
1486  {
1487  double mass_flux = 0;
1488  double water_flux = 0;
1489  for (unsigned int k=0; k<qsize; k++)
1490  water_flux += arma::dot(ad_coef[sbi][k], fe_values.normal_vector(k))*fe_values.JxW(k);
1491  water_flux /= side->measure();
1492 
1493  for (unsigned int k=0; k<qsize; k++)
1494  {
1495  double conc = 0;
1496  conc_grad.zeros();
1497  for (unsigned int i=0; i<ndofs; i++)
1498  {
1499  conc += fe_values.shape_value(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1500  conc_grad += fe_values.shape_grad(i,k)*ls[sbi]->get_solution_array()[dof_indices[i]-feo->dh()->loffset()];
1501  }
1502 
1503  // flux due to advection
1504  mass_flux += water_flux*conc*fe_values.JxW(k);
1505 
1506  if (bc_type[sbi] == EqData::dirichlet || (bc_type[sbi] == EqData::inflow && water_flux*side->measure() < -mh_dh->precision()))
1507  {
1508  // flux due to diffusion
1509  mass_flux -= arma::dot(dif_coef[sbi][k]*conc_grad,fe_values.normal_vector(k))*fe_values.JxW(k);
1510 
1511  // the penalty term has to be added otherwise the mass balance will not hold
1512  mass_flux -= gamma[sbi][side->cond_idx()]*(bc_values[k][sbi] - conc)*fe_values.JxW(k);
1513  }
1514 
1515  }
1516  bcd_balance[sbi][bc_region_idx] += mass_flux;
1517  if (mass_flux >= 0) bcd_plus_balance[sbi][bc_region_idx] += mass_flux;
1518  else bcd_minus_balance[sbi][bc_region_idx] += mass_flux;
1519  }
1520  }
1521 
1522 }
1523 
1524 template<class Model>
1526 {
1527  calc_elem_sources<1>(mass, src_balance);
1528  calc_elem_sources<2>(mass, src_balance);
1529  calc_elem_sources<3>(mass, src_balance);
1530 }
1531 
1532 template<class Model>
1533 template<unsigned int dim>
1535 {
1536  FEValues<dim,3> fe_values(*feo->mapping<dim>(), *feo->q<dim>(), *feo->fe<dim>(),
1538  const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim>()->size();
1539  unsigned int dof_indices[ndofs];
1540  vector<arma::vec> sources_conc(qsize), sources_density(qsize), sources_sigma(qsize);
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:237
#define DBGMSG(...)
Definition: global_defs.h:195
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:75
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:61
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:195
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:365
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:77
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).
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:78
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:120
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:180
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:95
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:104
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:227
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:274
unsigned int edge_idx_
Definition: boundaries.h:94
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:95
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:83
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:66
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:43
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:138
Definition: system.hh:75
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:313
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:449
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:83
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:72
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:70
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
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()
~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:228
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