45 using namespace Input::Type;
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");
57 .
add_value(inflow,
"inflow",
"Dirichlet BC on inflow and homogeneous Neumann BC on outflow.")
64 Model::ModelEqData::get_output_selection_input_type(
"DG",
"DG solver")
65 .copy_values(EqData().make_output_field_selection(
"").close())
70 = Model::get_input_type(
"DG",
"DG solver")
72 "Linear solver for MH problem.")
75 "Variant of interior penalty discontinuous Galerkin method.")
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.");
120 xprintf(
PrgErr,
"Unsupported polynomial order %d for finite elements in TransportDG ", fe_order);
140 dh_->distribute_dofs(*fe1_, *fe2_, *fe3_);
186 template<
class Model>
190 .
name(
"fracture_sigma")
192 "Coefficient of diffusive transfer through fractures (for each substance).")
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.")
207 "Boundary condition type, possible values: inflow, dirichlet, neumann, robin.")
215 .
description(
"Flux in Neumann boundary condition.")
220 .
name(
"bc_robin_sigma")
221 .
description(
"Conductivity coefficient in Robin boundary condition.")
230 template<
class Model>
237 static_assert(std::is_base_of<AdvectionDiffusionModel, Model>::value,
"");
253 data_.set_mesh(init_mesh);
265 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
273 int qsize = max(
feo->
q<0>()->size(), max(
feo->
q<1>()->size(), max(
feo->
q<2>()->size(),
feo->
q<3>()->size())));
278 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
285 for (
int sd=0; sd<max_edg_sides; sd++)
289 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
300 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
310 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
328 for (
unsigned int sbi = 0; sbi <
n_subst_; sbi++) {
338 template<
class Model>
344 if (feo->dh()->el_ds()->myp() == 0)
346 for (
unsigned int i=0; i<
n_subst_; i++)
348 VecDestroy(&output_vec[i]);
349 delete[] output_solution[i];
353 for (
unsigned int i=0; i<
n_subst_; i++)
356 MatDestroy(&stiffness_matrix[i]);
360 delete[] stiffness_matrix;
367 delete output_stream;
371 template<
class Model>
375 VecScatter output_scatter;
377 for (
unsigned int sbi=0; sbi<
n_subst_; 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));
391 template<
class Model>
394 data_.set_time(*
time_);
397 set_initial_condition();
398 for (
unsigned int sbi = 0; sbi <
n_subst_; sbi++)
399 ( (
LinSys_PETSC *)ls[sbi] )->set_initial_guess_nonzero();
405 template<
class Model>
414 data_.set_time(*
time_);
418 if (!allocation_done)
421 ls_dt->start_allocation();
422 assemble_mass_matrix();
426 for (
unsigned int i=0; i<
n_subst_; i++)
428 ls[i]->start_allocation();
429 stiffness_matrix[i] = NULL;
432 assemble_stiffness_matrix();
434 set_boundary_conditions();
436 allocation_done =
true;
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());
450 if (stiffness_matrix[0] == NULL
452 || Model::flux_changed)
456 for (
unsigned int i=0; i<
n_subst_; i++)
458 ls[i]->start_add_assembly();
459 ls[i]->mat_zero_entries();
461 assemble_stiffness_matrix();
462 for (
unsigned int i=0; i<
n_subst_; i++)
464 ls[i]->finish_assembly();
466 if (stiffness_matrix[i] == NULL)
467 MatConvert(*( ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &stiffness_matrix[i]);
469 MatCopy(*( ls[i]->get_matrix() ), stiffness_matrix[i], DIFFERENT_NONZERO_PATTERN);
476 || Model::flux_changed)
478 for (
unsigned int i=0; i<
n_subst_; i++)
480 ls[i]->start_add_assembly();
481 ls[i]->rhs_zero_entries();
484 set_boundary_conditions();
485 for (
unsigned int i=0; i<
n_subst_; i++)
487 ls[i]->finish_assembly();
489 VecDuplicate(*( ls[i]->get_rhs() ), &rhs[i]);
490 VecCopy(*( ls[i]->get_rhs() ), rhs[i]);
494 Model::flux_changed =
false;
515 for (
unsigned int i=0; i<
n_subst_; i++)
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);
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]);
542 template<
class Model>
549 Model::flux_changed =
true;
554 template<
class Model>
562 output_vector_gather();
564 data_.output(output_stream);
565 output_stream->write_time_frame();
574 template<
class Model>
578 assemble_mass_matrix<1>();
579 assemble_mass_matrix<2>();
580 assemble_mass_matrix<3>();
585 template<
class Model>
template<
unsigned int dim>
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];
594 for (
unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
597 if (cell->dim() != dim)
continue;
599 fe_values.reinit(cell);
600 feo->dh()->get_dof_indices(cell, dof_indices);
603 Model::compute_mass_matrix_coefficient(fe_values.point_list(), ele_acc, mm_coef);
606 for (
unsigned int i=0; i<ndofs; i++)
608 for (
unsigned int j=0; j<ndofs; j++)
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);
616 ls_dt->mat_set_values(ndofs, (
int *)dof_indices, ndofs, (
int *)dof_indices, local_mass_matrix);
623 template<
class Model>
628 assemble_volume_integrals<1>();
629 assemble_volume_integrals<2>();
630 assemble_volume_integrals<3>();
634 assemble_fluxes_boundary<1>();
635 assemble_fluxes_boundary<2>();
636 assemble_fluxes_boundary<3>();
640 assemble_fluxes_element_element<1>();
641 assemble_fluxes_element_element<2>();
642 assemble_fluxes_element_element<3>();
646 assemble_fluxes_element_side<1>();
647 assemble_fluxes_element_side<2>();
648 assemble_fluxes_element_side<3>();
655 template<
class Model>
656 template<
unsigned int dim>
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];
667 PetscScalar local_matrix[ndofs*ndofs];
670 for (
unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
673 if (cell->dim() != dim)
continue;
678 feo->dh()->get_dof_indices(cell, dof_indices);
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);
685 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
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;
691 for (
unsigned int k=0; k<qsize; k++)
693 for (
unsigned int i=0; i<ndofs; i++)
696 double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.
shape_grad(i,k));
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))
704 ls[sbi]->mat_set_values(ndofs, (
int *)dof_indices, ndofs, (
int *)dof_indices, local_matrix);
710 template<
class Model>
720 template<
class Model>
721 template<
unsigned int dim>
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];
733 for (
unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
736 if (cell->dim() != dim)
continue;
738 fe_values.reinit(cell);
739 feo->dh()->get_dof_indices(cell, dof_indices);
741 Model::compute_source_coefficients(fe_values.point_list(), cell->element_accessor(), sources_conc, sources_density, sources_sigma);
744 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
746 for (
unsigned int i=0; i<ndofs; i++)
750 for (
unsigned int k=0; k<qsize; k++)
752 source = (sources_density[k][sbi] + sources_conc[k][sbi]*sources_sigma[k][sbi])*fe_values.JxW(k);
754 for (
unsigned int i=0; i<ndofs; i++)
755 local_rhs[i] += source*fe_values.shape_value(i,k);
757 ls[sbi]->rhs_set_values(ndofs, (
int *)dof_indices, local_rhs);
764 template<
class Model>
765 template<
unsigned int dim>
771 const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size(),
772 n_max_sides = ad_coef_edg.size();
774 PetscScalar local_matrix[ndofs*ndofs];
777 double gamma_l, omega[2], transport_flux;
779 for (
unsigned int sid=0; sid<n_max_sides; sid++)
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>(),
787 for (
unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
792 for (
int sid=0; sid<edg->
n_sides; sid++)
796 feo->dh()->get_dof_indices(cell, side_dof_indices[sid]);
797 fe_values[sid]->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);
805 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
808 for (
int sid=0; sid<edg->
n_sides; sid++)
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);
816 for (
int s1=0; s1<edg->
n_sides; s1++)
818 for (
int s2=s1+1; s2<edg->
n_sides; s2++)
824 arma::vec3 nv = fe_values[s1]->normal_vector(0);
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);
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))
838 for (
int n=0; n<2; n++)
842 for (
int m=0; m<2; m++)
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;
848 for (
unsigned int k=0; k<qsize; k++)
850 double flux_times_JxW = transport_flux*fe_values[0]->JxW(k);
851 double gamma_times_JxW = gamma_l*fe_values[0]->JxW(k);
853 for (
unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
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;
860 for (
unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
862 int index = i*fe_values[sd[m]]->n_dofs()+j;
865 local_matrix[index] += flux_JxW_jump_i*
AVERAGE(j,k,m);
868 local_matrix[index] += gamma_JxW_jump_i*
JUMP(j,k,m);
871 local_matrix[index] -=
WAVERAGE(j,k,m)*JxW_jump_i;
872 local_matrix[index] -=
JUMP(j,k,m)*JxW_var_wavg_i;
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);
887 for (
unsigned int i=0; i<n_max_sides; i++)
890 delete[] side_dof_indices[i];
895 template<
class Model>
896 template<
unsigned int dim>
899 FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
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];
908 arma::vec dg_penalty;
912 for (
unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
915 if (edg->
n_sides > 1)
continue;
917 if (edg->
side(0)->
dim() != dim-1)
continue;
919 if (edg->
side(0)->
cond() == NULL)
continue;
924 feo->dh()->get_dof_indices(cell, side_dof_indices);
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);
934 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
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;
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();
947 if ((bc_type[sbi] == EqData::dirichlet) || (bc_type[sbi] == EqData::inflow ))
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;
957 for (
unsigned int k=0; k<qsize; k++)
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);
963 flux_times_JxW = transport_flux*fe_values_side.
JxW(k);
965 for (
unsigned int i=0; i<ndofs; i++)
966 for (
unsigned int j=0; j<ndofs; j++)
969 local_matrix[i*ndofs+j] += flux_times_JxW*fe_values_side.
shape_value(i,k)*fe_values_side.
shape_value(j,k);
972 if (bc_type[sbi] == EqData::dirichlet || (bc_type[sbi] == EqData::inflow && side_flux < -mh_dh->precision()))
975 )*fe_values_side.
JxW(k);
978 ls[sbi]->mat_set_values(ndofs, (
int *)side_dof_indices, ndofs, (
int *)side_dof_indices, local_matrix);
984 template<
class Model>
985 template<
unsigned int dim>
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>(),
996 FEValues<dim-1,3> fv_rt(*feo->mapping<dim-1>(), *feo->q<dim-1>(), *feo->fe_rt<dim-1>(),
1000 const unsigned int ndofs = feo->fe<dim>()->n_dofs();
1001 const unsigned int qsize = feo->q<dim-1>()->size();
1002 unsigned int side_dof_indices[2*ndofs], n_dofs[2];
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];
1011 fv_sb[0] = &fe_values_vb;
1012 fv_sb[1] = &fe_values_side;
1015 for (
unsigned int inb=0; inb<feo->dh()->n_loc_nb(); inb++)
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();
1027 feo->dh()->get_dof_indices(cell, side_dof_indices+n_dofs[0]);
1029 n_dofs[1] = fv_sb[1]->n_dofs();
1033 element_id[0] = cell_sub.
index();
1034 element_id[1] = cell.
index();
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);
1048 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
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;
1055 for (
unsigned int k=0; k<qsize; k++)
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]);
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]);
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);
1077 for (
int n=0; n<2; n++)
1079 if (!feo->dh()->el_is_local(element_id[n]))
continue;
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);
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);
1099 template<
class Model>
1103 set_boundary_conditions<1>();
1104 set_boundary_conditions<2>();
1105 set_boundary_conditions<3>();
1110 template<
class Model>
1111 template<
unsigned int dim>
1114 FESideValues<dim,3> fe_values_side(*feo->mapping<dim>(), *feo->q<dim-1>(), *feo->fe<dim>(),
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];
1124 for (
unsigned int i=0; i<qsize; i++)
1126 bc_values[i].resize(n_subst_);
1127 bc_fluxes[i].resize(n_subst_);
1128 bc_sigma[i].resize(n_subst_);
1131 for (
unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1134 if (edg->
n_sides > 1)
continue;
1135 if (edg->
side(0)->
dim() != dim-1)
continue;
1137 if (edg->
side(0)->
cond() == NULL)
continue;
1143 arma::uvec bc_type = data_.bc_type.value(side->
cond()->
element()->
centre(), ele_acc);
1145 fe_values_side.reinit(cell, side->
el_idx());
1146 fsv_rt.reinit(cell, side->
el_idx());
1147 calculate_velocity(cell, velocity, fsv_rt);
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);
1154 feo->dh()->get_dof_indices(cell, side_dof_indices);
1156 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
1158 for (
unsigned int i=0; i<ndofs; i++) local_rhs[i] = 0;
1160 for (
unsigned int k=0; k<qsize; k++)
1166 || (bc_type[sbi] == EqData::dirichlet))
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));
1171 else if (bc_type[sbi] == EqData::neumann)
1173 bc_term = bc_fluxes[k][sbi]*fe_values_side.JxW(k);
1175 else if (bc_type[sbi] == EqData::robin)
1177 bc_term = bc_sigma[k][sbi]*bc_values[k][sbi]*fe_values_side.JxW(k);
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));
1184 ls[sbi]->rhs_set_values(ndofs, (
int *)side_dof_indices, local_rhs);
1196 template<
class Model>
1197 template<
unsigned int dim>
1201 for (
unsigned int i=0; i<cell->n_nodes(); i++)
1202 node_nums[cell->node[i]] = i;
1206 for (
unsigned int k=0; k<fv.
n_points(); k++)
1208 velocity[k].zeros();
1209 for (
unsigned int sid=0; sid<cell->n_sides(); sid++)
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)];
1224 template<
class Model>
1233 const double alpha1,
1234 const double alpha2,
1237 double &transport_flux)
1241 double local_alpha = 0;
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)));
1259 double pflux = 0, nflux = 0;
1260 for (
int i=0; i<edg.
n_sides; i++)
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);
1274 transport_flux = fluxes[s1];
1278 gamma = 0.5*fabs(transport_flux);
1282 local_alpha = max(alpha1, alpha2);
1290 for (
int k=0; k<K_size; k++)
1291 delta[0] += dot(K1[k]*normal_vector,normal_vector);
1294 gamma += local_alpha/h*(delta[0]);
1300 for (
int k=0; k<K_size; k++)
1302 delta[0] += dot(K1[k]*normal_vector,normal_vector);
1303 delta[1] += dot(K2[k]*normal_vector,normal_vector);
1308 double delta_sum = delta[0] + delta[1];
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);
1317 for (
int i=0; i<2; i++) omega[i] = 0;
1326 template<
class Model>
1335 double delta = 0, h = 0;
1338 if (side->
dim() == 0)
1344 for (
unsigned int i=0; i<side->
n_nodes(); i++)
1345 for (
unsigned int j=i+1; j<side->
n_nodes(); j++)
1350 for (
int k=0; k<K_size; k++)
1351 delta += dot(K[k]*normal_vector,normal_vector);
1354 gamma = 0.5*fabs(flux) + alpha/h*delta;
1361 template<
class Model>
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>();
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>();
1377 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
1379 ls[sbi]->finish_assembly();
1385 template<
class Model>
1386 template<
unsigned int dim>
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];
1396 for (
unsigned int k=0; k<qsize; k++)
1397 init_values[k].resize(n_subst_);
1399 for (
unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
1402 if (elem->dim() != dim)
continue;
1405 feo->dh()->get_dof_indices(elem, dof_indices);
1408 Model::compute_init_cond(fe_values.
point_list(), ele_acc, init_values);
1410 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
1412 for (
unsigned int i=0; i<ndofs; i++)
1415 for (
unsigned int j=0; j<ndofs; j++)
1416 matrix[i*ndofs+j] = 0;
1419 for (
unsigned int k=0; k<qsize; k++)
1421 double rhs_term = init_values[k](sbi)*fe_values.
JxW(k);
1423 for (
unsigned int i=0; i<ndofs; i++)
1425 for (
unsigned int j=0; j<ndofs; j++)
1431 ls[sbi]->set_values(ndofs, (
int *)dof_indices, ndofs, (
int *)dof_indices, matrix, rhs);
1437 template<
class Model>
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);
1445 template<
class Model>
1446 template<
unsigned int dim>
1453 const unsigned int ndofs = feo->fe<dim>()->n_dofs(), qsize = feo->q<dim-1>()->size();
1454 unsigned int dof_indices[ndofs];
1459 for (
unsigned int k=0; k<qsize; k++)
1460 bc_values[k].resize(n_subst_);
1462 for (
unsigned int iedg=0; iedg<feo->dh()->n_loc_edges(); iedg++)
1465 if (edg->
n_sides > 1)
continue;
1466 if (edg->
side(0)->
dim() != dim-1)
continue;
1468 if (edg->
side(0)->
cond() == NULL)
continue;
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);
1485 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
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();
1493 for (
unsigned int k=0; k<qsize; k++)
1497 for (
unsigned int i=0; i<ndofs; i++)
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()];
1504 mass_flux += water_flux*conc*fe_values.JxW(k);
1506 if (bc_type[sbi] == EqData::dirichlet || (bc_type[sbi] == EqData::inflow && water_flux*side->
measure() < -
mh_dh->
precision()))
1509 mass_flux -= arma::dot(dif_coef[sbi][k]*conc_grad,fe_values.normal_vector(k))*fe_values.JxW(k);
1512 mass_flux -= gamma[sbi][side->
cond_idx()]*(bc_values[k][sbi] - conc)*fe_values.JxW(k);
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;
1524 template<
class Model>
1527 calc_elem_sources<1>(mass, src_balance);
1528 calc_elem_sources<2>(mass, src_balance);
1529 calc_elem_sources<3>(mass, src_balance);
1532 template<
class Model>
1533 template<
unsigned int dim>
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;
1543 for (
unsigned int i_cell=0; i_cell<feo->dh()->el_ds()->lsize(); i_cell++)
1546 if (elem->dim() != dim)
continue;
1548 Region r = elem->element_accessor().region();
1550 unsigned int region_idx = r.
bulk_idx();
1553 fe_values.reinit(elem);
1554 feo->dh()->get_dof_indices(elem, dof_indices);
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);
1559 for (
unsigned int sbi=0; sbi<
n_subst_; sbi++)
1564 for (
unsigned int k=0; k<qsize; k++)
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()];
1570 mass_sum += mm_coef[k]*conc*fe_values.JxW(k);
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);
1576 mass[sbi][region_idx] += mass_sum;
1577 src_balance[sbi][region_idx] += sources_sum;
static Input::Type::Selection output_selection
void calculate(double time)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
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.
void assemble_fluxes_element_side()
Assembles the fluxes between elements of different dimensions.
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
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)
Transport with dispersion implemented using discontinuous Galerkin method.
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
double distance(const Node &n2) const
void assemble_mass_matrix()
Assembles the mass matrix.
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...
const arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
TimeMark::Type type_output()
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
const MH_DofHandler * mh_dh
MassBalance * mass_balance_
object for calculation and writing the mass balance to file.
void set_initial_condition()
Sets the initial condition.
Class FEValues calculates finite element data on the actual cells such as shape function values...
bool is_current(const TimeMark::Type &mask) const
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
static Input::Type::Record input_type
const RegionDB & region_db() const
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.
vector< Boundary > boundary_
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...
static TimeMarks & marks()
MassBalance * mass_balance()
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.
unsigned int boundary_idx() const
Returns index of the region in the boundary set.
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.
Transformed quadrature points.
void output_data()
Postprocesses the solution and writes to output file.
static OutputTime * create_output_stream(const Input::Record &in_rec)
This method write all registered data to output streams.
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
FEObjects * feo
Finite element objects.
static constexpr Mask in_time_term
A field is part of time term of the equation.
FieldCommon & input_default(const string &input_default)
Raviart-Thomas element of order 0.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
FieldCommon & input_selection(const Input::Type::Selection *element_selection)
Mat mass_matrix
The mass matrix.
#define START_TIMER(tag)
Starts a timer with specified tag.
FiniteElement< dim, spacedim > * fe
The used finite element.
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Mat * stiffness_matrix
The stiffness matrix.
void output(double time)
Write computed fields to file.
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.
unsigned int n_subst_
Number of transported substances.
Affine mapping between reference and actual cell.
Shape function gradients.
void mark_output_times(const TimeGovernor &tg)
Vec * rhs
Vector of right hand side.
void set_solution(double *sol_array)
BCField< 3, FieldValue< 3 >::Vector > bc_flux
Flux in Neumann or Robin b.c.
static Input::Type::Selection bc_type_selection
DOFHandlerMultiDim * dh()
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.
std::vector< string > subst_names_
Names of transported substances.
void output_data() override
Write computed fields.
vector< Neighbour > vb_neighbours_
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...
FieldCommon & description(const string &description)
EqData data_
Field data for model parameters.
bool allocation_done
Indicates whether matrices have been preallocated.
unsigned int cond_idx() const
std::vector< std::vector< double > > gamma
Penalty parameters.
unsigned int max_edge_sides(unsigned int dim) const
Discontinuous Galerkin method for equation of transport with dispersion.
ElementFullIter element() const
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.
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.
arma::vec3 centre() const
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
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...
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.)
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
FullIter full_iter(Iter it)
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Abstract linear system class.
Definitions of particular quadrature rules on simplices.
FieldCommon & name(const string &name)
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.
#define END_TIMER(tag)
Ends a timer with specified tag.
Discontinuous Galerkin method for equation of transport with dispersion.
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
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.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
const Node * node(unsigned int i) const
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
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
Field< 3, FieldValue< 3 >::Vector > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
#define JUMP(i, k, side_id)
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.
BCField< 3, FieldValue< 3 >::Vector > bc_robin_sigma
Transition coefficient in Robin b.c.
ElementAccessor< 3 > element_accessor()
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)
ElementVector element
Vector of elements of the mesh.
unsigned int n_nodes() const
Transformed quadrature weights.
Calculates finite element data on a side.