54 return Selection(
"DG_variant",
"Type of penalty term.")
55 .
add_value(non_symmetric,
"non-symmetric",
"non-symmetric weighted interior penalty DG method")
56 .
add_value(incomplete,
"incomplete",
"incomplete weighted interior penalty DG method")
57 .
add_value(symmetric,
"symmetric",
"symmetric weighted interior penalty DG method")
80 std::string equation_name = std::string(Model::ModelEqData::name()) +
"_DG";
81 return Model::get_input_type(
"DG",
"Discontinuous Galerkin (DG) solver")
83 "Solver for the linear system.")
86 .make_field_descriptor_type(equation_name)),
88 "Input fields of the equation.")
90 "Variant of the interior penalty discontinuous Galerkin method.")
92 "Polynomial order for the finite element in DG method (order 0 is suitable if there is no diffusion/dispersion).")
94 EqData().output_fields.make_output_type(equation_name,
""),
95 IT::Default(
"{ \"fields\": [ " + Model::ModelEqData::default_output_field() +
"] }"),
96 "Specification of output fields and output times.")
100 template<
class Model>
102 Input::register_class< TransportDG<Model>,
Mesh &,
const Input::Record>(std::string(Model::ModelEqData::name()) +
"_DG") +
110 unsigned int q_order;
112 q_order = 2*fe_order;
131 ds_ = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, fe0_, fe1_, fe2_, fe3_);
132 dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
134 dh_->distribute_dofs(ds_);
178 template<
class Model>
182 .
name(
"fracture_sigma")
184 "Coefficient of diffusive transfer through fractures (for each substance).")
192 "Penalty parameter influencing the discontinuity of the solution (for each substance). " 193 "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
206 .
description(
"Subdomain ids of the domain decomposition.");
214 template<
class Model>
216 : Model(init_mesh, in_rec),
230 data_.set_mesh(init_mesh);
239 Model::init_from_input(in_rec);
248 template<
class Model>
251 data_.set_components(Model::substances_.names());
255 gamma.resize(Model::n_substances());
256 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
257 gamma[sbi].resize(Model::mesh_->boundary_.size());
260 int qsize = max(
feo->
q<0>()->size(), max(
feo->
q<1>()->size(), max(
feo->
q<2>()->size(),
feo->
q<3>()->size())));
261 int max_edg_sides = max(Model::mesh_->max_edge_sides(1), max(Model::mesh_->max_edge_sides(2), Model::mesh_->max_edge_sides(3)));
263 ret_coef.resize(Model::n_substances());
266 ad_coef.resize(Model::n_substances());
267 dif_coef.resize(Model::n_substances());
268 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
276 for (
int sd=0; sd<max_edg_sides; sd++)
280 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
288 data_.output_field.set_components(Model::substances_.names());
289 data_.output_field.set_mesh(*Model::mesh_);
292 data_.output_field.setup_components();
293 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
298 data_.output_field[sbi].set_field(Model::mesh_->region_db().get_region_set(
"ALL"), output_field_ptr, 0);
305 std::string petsc_default_opts;
306 if (
feo->
dh()->distr()->np() == 1)
307 petsc_default_opts =
"-ksp_type bcgs -pc_type ilu -pc_factor_levels 2 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
309 petsc_default_opts =
"-ksp_type bcgs -ksp_diagonal_scale_fix -pc_type asm -pc_asm_overlap 4 -sub_pc_type ilu -sub_pc_factor_levels 3 -sub_pc_factor_fill 6.0";
312 ls =
new LinSys*[Model::n_substances()];
317 mass_matrix.resize(Model::n_substances(),
nullptr);
318 rhs.resize(Model::n_substances(),
nullptr);
319 mass_vec.resize(Model::n_substances(),
nullptr);
320 ret_vec.resize(Model::n_substances(),
nullptr);
322 for (
unsigned int sbi = 0; sbi < Model::n_substances(); sbi++) {
329 solution_elem_[sbi] =
new double[Model::mesh_->get_el_ds()->lsize()];
336 Model::balance_->allocate(
feo->
dh()->distr()->lsize(),
337 max(
feo->
fe<1>()->n_dofs(), max(
feo->
fe<2>()->n_dofs(),
feo->
fe<3>()->n_dofs())));
342 template<
class Model>
347 if (
gamma.size() > 0) {
350 for (
unsigned int i=0; i<Model::n_substances(); i++)
382 template<
class Model>
386 data_.mark_input_times( *(Model::time_) );
388 std::stringstream ss;
396 for (
unsigned int sbi = 0; sbi < Model::n_substances(); sbi++)
403 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
405 Model::balance_->calculate_instant(Model::subst_idx[sbi],
ls[sbi]->
get_solution());
415 template<
class Model>
419 for (
unsigned int i=0; i<Model::n_substances(); i++)
435 for (
unsigned int i=0; i<Model::n_substances(); i++)
446 template<
class Model>
451 Model::time_->next_time();
452 Model::time_->view(
"TDG");
461 for (
unsigned int i=0; i<Model::n_substances(); i++)
468 for (
unsigned int i=0; i<Model::n_substances(); i++)
478 MatConvert(*(
ls_dt[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &
mass_matrix[i]);
481 MatCopy(*(
ls_dt[i]->get_matrix() ),
mass_matrix[i], DIFFERENT_NONZERO_PATTERN);
488 || Model::flux_changed)
492 for (
unsigned int i=0; i<Model::n_substances(); i++)
498 for (
unsigned int i=0; i<Model::n_substances(); i++)
503 MatConvert(*(
ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &
stiffness_matrix[i]);
512 || Model::flux_changed)
514 for (
unsigned int i=0; i<Model::n_substances(); i++)
521 for (
unsigned int i=0; i<Model::n_substances(); i++)
525 if (
rhs[i] ==
nullptr) VecDuplicate(*(
ls[i]->get_rhs() ), &
rhs[i]);
526 VecCopy(*(
ls[i]->get_rhs() ),
rhs[i]);
530 Model::flux_changed =
false;
551 for (
unsigned int i=0; i<Model::n_substances(); i++)
554 MatAXPY(m, 1./Model::time_->dt(),
mass_matrix[i], SUBSET_NONZERO_PATTERN);
557 VecDuplicate(
rhs[i], &w);
558 VecWAXPY(w, 1./Model::time_->dt(),
mass_vec[i],
rhs[i]);
577 template<
class Model>
581 unsigned int i_cell=0;
582 for (
auto cell :
feo->
dh()->own_range() )
589 n_dofs =
feo->
fe<1>()->n_dofs();
592 n_dofs =
feo->
fe<2>()->n_dofs();
595 n_dofs =
feo->
fe<3>()->n_dofs();
600 cell.get_dof_indices(dof_indices);
602 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
606 for (
unsigned int j=0; j<n_dofs; ++j)
607 solution_elem_[sbi][i_cell] +=
ls[sbi]->get_solution_array()[dof_indices[j]-
feo->
dh()->distr()->begin()];
618 template<
class Model>
631 Model::output_data();
634 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
635 Model::balance_->calculate_instant(Model::subst_idx[sbi],
ls[sbi]->
get_solution());
636 Model::balance_->output();
643 template<
class Model>
646 if (Model::balance_->cumulative())
648 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
650 Model::balance_->calculate_cumulative(Model::subst_idx[sbi],
ls[sbi]->
get_solution());
662 template<
class Model>
666 Model::balance_->start_mass_assembly(Model::subst_idx);
667 assemble_mass_matrix<1>();
668 assemble_mass_matrix<2>();
669 assemble_mass_matrix<3>();
670 Model::balance_->finish_mass_assembly(Model::subst_idx);
675 template<
class Model>
template<
unsigned int dim>
679 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
681 PetscScalar local_mass_matrix[ndofs*ndofs], local_retardation_balance_vector[ndofs];
685 for (
auto cell :
feo->
dh()->own_range() )
687 if (cell.dim() != dim)
continue;
691 cell.get_dof_indices(dof_indices);
696 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
699 for (
unsigned int i=0; i<ndofs; i++)
701 for (
unsigned int j=0; j<ndofs; j++)
703 local_mass_matrix[i*ndofs+j] = 0;
704 for (
unsigned int k=0; k<qsize; k++)
709 for (
unsigned int i=0; i<ndofs; i++)
711 local_mass_balance_vector[i] = 0;
712 local_retardation_balance_vector[i] = 0;
713 for (
unsigned int k=0; k<qsize; k++)
720 Model::balance_->add_mass_matrix_values(Model::subst_idx[sbi], elm.
region().
bulk_idx(), dof_indices, local_mass_balance_vector);
721 ls_dt[sbi]->
mat_set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), local_mass_matrix);
722 VecSetValues(
ret_vec[sbi], ndofs, &(dof_indices[0]), local_retardation_balance_vector, ADD_VALUES);
730 template<
class Model>
735 assemble_volume_integrals<1>();
736 assemble_volume_integrals<2>();
737 assemble_volume_integrals<3>();
741 assemble_fluxes_boundary<1>();
742 assemble_fluxes_boundary<2>();
743 assemble_fluxes_boundary<3>();
747 assemble_fluxes_element_element<1>();
748 assemble_fluxes_element_element<2>();
749 assemble_fluxes_element_element<3>();
753 assemble_fluxes_element_side<1>();
754 assemble_fluxes_element_side<2>();
755 assemble_fluxes_element_side<3>();
762 template<
class Model>
763 template<
unsigned int dim>
770 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
774 PetscScalar local_matrix[ndofs*ndofs];
777 for (
auto cell :
feo->
dh()->local_range() )
779 if (!cell.is_own())
continue;
780 if (cell.dim() != dim)
continue;
785 cell.get_dof_indices(dof_indices);
789 Model::compute_sources_sigma(fe_values.
point_list(), elm, sources_sigma);
792 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
794 for (
unsigned int i=0; i<ndofs; i++)
795 for (
unsigned int j=0; j<ndofs; j++)
796 local_matrix[i*ndofs+j] = 0;
798 for (
unsigned int k=0; k<qsize; k++)
800 for (
unsigned int i=0; i<ndofs; i++)
803 double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.
shape_grad(i,k));
805 for (
unsigned int j=0; j<ndofs; j++)
806 local_matrix[i*ndofs+j] += (arma::dot(Kt_grad_i, fe_values.
shape_grad(j,k))
811 ls[sbi]->
mat_set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), local_matrix);
817 template<
class Model>
821 Model::balance_->start_source_assembly(Model::subst_idx);
825 Model::balance_->finish_source_assembly(Model::subst_idx);
829 template<
class Model>
830 template<
unsigned int dim>
835 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
841 PetscScalar local_rhs[ndofs];
846 for (
auto cell :
feo->
dh()->own_range() )
848 if (cell.dim() != dim)
continue;
852 cell.get_dof_indices(dof_indices);
853 cell.get_loc_dof_indices(loc_dof_indices);
855 Model::compute_source_coefficients(fe_values.
point_list(), elm, sources_conc, sources_density, sources_sigma);
858 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
860 fill_n(local_rhs, ndofs, 0);
861 local_source_balance_vector.assign(ndofs, 0);
862 local_source_balance_rhs.assign(ndofs, 0);
865 for (
unsigned int k=0; k<qsize; k++)
867 source = (sources_density[sbi][k] + sources_conc[sbi][k]*sources_sigma[sbi][k])*fe_values.
JxW(k);
869 for (
unsigned int i=0; i<ndofs; i++)
874 for (
unsigned int i=0; i<ndofs; i++)
876 for (
unsigned int k=0; k<qsize; k++)
877 local_source_balance_vector[i] -= sources_sigma[sbi][k]*fe_values.
shape_value(i,k)*fe_values.
JxW(k);
879 local_source_balance_rhs[i] += local_rhs[i];
881 Model::balance_->add_source_values(Model::subst_idx[sbi], elm.
region().
bulk_idx(), loc_dof_indices,
882 local_source_balance_vector, local_source_balance_rhs);
892 double h_max = 0, h_min = numeric_limits<double>::infinity();
893 for (
unsigned int i=0; i<e->
n_nodes(); i++)
894 for (
unsigned int j=i+1; j<e->
n_nodes(); j++)
904 template<
class Model>
905 template<
unsigned int dim>
911 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim-1>()->size(),
914 PetscScalar local_matrix[ndofs*ndofs];
917 double gamma_l, omega[2], transport_flux, delta[2], delta_sum;
918 double aniso1, aniso2;
920 for (
unsigned int sid=0; sid<n_max_sides; sid++)
930 if (dh_cell.dim() != dim)
continue;
931 for(
DHCellSide cell_side : dh_cell.side_range() )
933 if (cell_side.n_edge_sides() < 2)
continue;
934 bool unique_edge = (cell_side.edge_sides().begin()->element().idx() != dh_cell.elm_idx());
935 if ( unique_edge )
continue;
939 auto dh_edge_cell =
feo->
dh()->cell_accessor_from_element( edge_side.elem_idx() );
941 dh_edge_cell.get_dof_indices(side_dof_indices[sid]);
942 fe_values[sid]->reinit(cell, edge_side.side_idx());
943 fsv_rt.
reinit(cell, edge_side.side_idx());
945 Model::compute_advection_diffusion_coefficients(fe_values[sid]->point_list(), side_velocity[sid], cell,
ad_coef_edg[sid],
dif_coef_edg[sid]);
946 dg_penalty[sid].resize(Model::n_substances());
947 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
951 arma::vec3 normal_vector = fe_values[0]->normal_vector(0);
954 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
957 double pflux = 0, nflux = 0;
962 for (
unsigned int k=0; k<qsize; k++)
963 fluxes[sid] += arma::dot(
ad_coef_edg[sid][sbi][k], fe_values[sid]->normal_vector(k))*fe_values[sid]->JxW(k);
964 fluxes[sid] /= edge_side.measure();
966 pflux += fluxes[sid];
968 nflux += fluxes[sid];
979 if (s2<=s1)
continue;
980 ASSERT(edge_side1.is_valid()).error(
"Invalid side of edge.");
982 arma::vec3 nv = fe_values[s1]->normal_vector(0);
986 if (fluxes[s2] > 0 && fluxes[s1] < 0)
987 transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
988 else if (fluxes[s2] < 0 && fluxes[s1] > 0)
989 transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
993 gamma_l = 0.5*fabs(transport_flux);
997 for (
unsigned int k=0; k<qsize; k++)
999 delta[0] += dot(
dif_coef_edg[s1][sbi][k]*normal_vector,normal_vector);
1000 delta[1] += dot(
dif_coef_edg[s2][sbi][k]*normal_vector,normal_vector);
1005 delta_sum = delta[0] + delta[1];
1008 if (fabs(delta_sum) > 0)
1010 omega[0] = delta[1]/delta_sum;
1011 omega[1] = delta[0]/delta_sum;
1012 double local_alpha = max(dg_penalty[s1][sbi], dg_penalty[s2][sbi]);
1013 double h = edge_side1.diameter();
1016 gamma_l += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
1019 for (
int i=0; i<2; i++) omega[i] = 0;
1022 int sd[2];
bool is_side_own[2];
1023 sd[0] = s1; is_side_own[0] = edge_side1.cell().is_own();
1024 sd[1] = s2; is_side_own[1] = edge_side2.cell().is_own();
1026 #define AVERAGE(i,k,side_id) (fe_values[sd[side_id]]->shape_value(i,k)*0.5) 1027 #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]) 1028 #define JUMP(i,k,side_id) ((side_id==0?1:-1)*fe_values[sd[side_id]]->shape_value(i,k)) 1031 for (
int n=0; n<2; n++)
1033 if (!is_side_own[n])
continue;
1035 for (
int m=0; m<2; m++)
1037 for (
unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
1038 for (
unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
1039 local_matrix[i*fe_values[sd[m]]->n_dofs()+j] = 0;
1041 for (
unsigned int k=0; k<qsize; k++)
1043 double flux_times_JxW = transport_flux*fe_values[0]->JxW(k);
1044 double gamma_times_JxW = gamma_l*fe_values[0]->JxW(k);
1046 for (
unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
1048 double flux_JxW_jump_i = flux_times_JxW*
JUMP(i,k,n);
1049 double gamma_JxW_jump_i = gamma_times_JxW*
JUMP(i,k,n);
1050 double JxW_jump_i = fe_values[0]->JxW(k)*
JUMP(i,k,n);
1053 for (
unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
1055 int index = i*fe_values[sd[m]]->n_dofs()+j;
1058 local_matrix[index] += flux_JxW_jump_i*
AVERAGE(j,k,m);
1061 local_matrix[index] += gamma_JxW_jump_i*
JUMP(j,k,m);
1064 local_matrix[index] -=
WAVERAGE(j,k,m)*JxW_jump_i;
1065 local_matrix[index] -=
JUMP(j,k,m)*JxW_var_wavg_i;
1069 ls[sbi]->
mat_set_values(fe_values[sd[n]]->n_dofs(), &(side_dof_indices[sd[n]][0]), fe_values[sd[m]]->n_dofs(), &(side_dof_indices[sd[m]][0]), local_matrix);
1082 for (
unsigned int i=0; i<n_max_sides; i++)
1084 delete fe_values[i];
1089 template<
class Model>
1090 template<
unsigned int dim>
1097 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim-1>()->size();
1099 PetscScalar local_matrix[ndofs*ndofs];
1106 for (
auto cell :
feo->
dh()->local_range() )
1108 if (!cell.is_own())
continue;
1109 for(
DHCellSide cell_side : cell.side_range() )
1111 Side side = cell_side.side();
1114 if (side.
dim() != dim-1)
continue;
1116 if (side.
cond() == NULL)
continue;
1119 cell.get_dof_indices(side_dof_indices);
1124 Model::compute_advection_diffusion_coefficients(fe_values_side.
point_list(), side_velocity, elm_acc,
ad_coef,
dif_coef);
1127 data_.cross_section.value_list(fe_values_side.
point_list(), elm_acc, csection);
1129 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1131 for (
unsigned int i=0; i<ndofs; i++)
1132 for (
unsigned int j=0; j<ndofs; j++)
1133 local_matrix[i*ndofs+j] = 0;
1137 double side_flux = 0;
1138 for (
unsigned int k=0; k<qsize; k++)
1139 side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.
normal_vector(k))*fe_values_side.
JxW(k);
1140 double transport_flux = side_flux/side.
measure();
1147 transport_flux += gamma_l;
1151 for (
unsigned int k=0; k<qsize; k++)
1153 double flux_times_JxW;
1157 flux_times_JxW = csection[k]*robin_sigma[k]*fe_values_side.
JxW(k);
1162 flux_times_JxW = (transport_flux + csection[k]*robin_sigma[k])*fe_values_side.
JxW(k);
1167 flux_times_JxW = transport_flux*fe_values_side.
JxW(k);
1169 for (
unsigned int i=0; i<ndofs; i++)
1171 for (
unsigned int j=0; j<ndofs; j++)
1174 local_matrix[i*ndofs+j] += flux_times_JxW*fe_values_side.
shape_value(i,k)*fe_values_side.
shape_value(j,k);
1180 )*fe_values_side.
JxW(k);
1185 ls[sbi]->
mat_set_values(ndofs, &(side_dof_indices[0]), ndofs, &(side_dof_indices[0]), local_matrix);
1192 template<
class Model>
1193 template<
unsigned int dim>
1197 if (dim == 1)
return;
1208 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs();
1209 const unsigned int qsize =
feo->
q<dim-1>()->size();
1210 int side_dof_indices[2*ndofs];
1212 unsigned int n_dofs[2], n_indices;
1216 PetscScalar local_matrix[4*ndofs*ndofs];
1217 double comm_flux[2][2];
1221 fv_sb[0] = &fe_values_vb;
1222 fv_sb[1] = &fe_values_side;
1226 for(
DHCellSide neighb_side : cell_lower_dim.neighb_sides() )
1229 if (cell_lower_dim.elm().dim() != dim-1)
continue;
1232 n_indices = cell_lower_dim.get_dof_indices(indices);
1233 for(
unsigned int i=0; i<n_indices; ++i) {
1234 side_dof_indices[i] = indices[i];
1236 fe_values_vb.
reinit(elm_lower_dim);
1237 n_dofs[0] = fv_sb[0]->n_dofs();
1239 DHCellAccessor cell_higher_dim =
feo->
dh()->cell_accessor_from_element( neighb_side.element().idx() );
1242 for(
unsigned int i=0; i<n_indices; ++i) {
1243 side_dof_indices[i+n_dofs[0]] = indices[i];
1245 fe_values_side.
reinit(elm_higher_dim, neighb_side.side_idx());
1246 n_dofs[1] = fv_sb[1]->n_dofs();
1249 bool own_element_id[2];
1250 own_element_id[0] = cell_lower_dim.is_own();
1251 own_element_id[1] = cell_higher_dim.
is_own();
1253 fsv_rt.
reinit(elm_higher_dim, neighb_side.side_idx());
1254 fv_rt.
reinit(elm_lower_dim);
1258 Model::compute_advection_diffusion_coefficients(fe_values_vb.
point_list(), velocity_higher, elm_higher_dim, ad_coef_edg[1],
dif_coef_edg[1]);
1259 data_.cross_section.value_list(fe_values_vb.
point_list(), elm_lower_dim, csection_lower);
1260 data_.cross_section.value_list(fe_values_vb.
point_list(), elm_higher_dim, csection_higher);
1262 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1264 for (
unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
1265 for (
unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
1266 local_matrix[i*(n_dofs[0]+n_dofs[1])+j] = 0;
1271 for (
unsigned int k=0; k<qsize; k++)
1282 2*csection_higher[k]*csection_higher[k]/(csection_lower[k]*csection_lower[k]);
1284 double transport_flux = arma::dot(ad_coef_edg[1][sbi][k], fe_values_side.
normal_vector(k));
1286 comm_flux[0][0] = (sigma-min(0.,transport_flux))*fv_sb[0]->JxW(k);
1287 comm_flux[0][1] = -(sigma-min(0.,transport_flux))*fv_sb[0]->JxW(k);
1288 comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1289 comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1291 for (
int n=0; n<2; n++)
1293 if (!own_element_id[n])
continue;
1295 for (
unsigned int i=0; i<n_dofs[n]; i++)
1296 for (
int m=0; m<2; m++)
1297 for (
unsigned int j=0; j<n_dofs[m]; j++)
1298 local_matrix[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
1299 comm_flux[m][n]*fv_sb[m]->shape_value(j,k)*fv_sb[n]->shape_value(i,k);
1302 ls[sbi]->
mat_set_values(n_dofs[0]+n_dofs[1], side_dof_indices, n_dofs[0]+n_dofs[1], side_dof_indices, local_matrix);
1313 template<
class Model>
1317 Model::balance_->start_flux_assembly(Model::subst_idx);
1318 set_boundary_conditions<1>();
1319 set_boundary_conditions<2>();
1320 set_boundary_conditions<3>();
1321 Model::balance_->finish_flux_assembly(Model::subst_idx);
1326 template<
class Model>
1327 template<
unsigned int dim>
1334 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim-1>()->size();
1336 double local_rhs[ndofs];
1338 PetscScalar local_flux_balance_rhs;
1342 bc_ref_values(qsize);
1346 for (
auto cell :
feo->
dh()->own_range() )
1348 if (cell.elm()->boundary_idx_ ==
nullptr)
continue;
1350 for (
unsigned int si=0; si<cell.elm()->n_sides(); si++)
1353 if (edg->
n_sides > 1)
continue;
1355 if (edg->
side(0)->
cond() == NULL)
continue;
1358 if (edg->
side(0)->
dim() != dim-1)
continue;
1366 Model::get_bc_type(ele_acc, bc_type);
1375 auto dh_cell =
feo->
dh()->cell_accessor_from_element( side.
element().
idx() );
1376 dh_cell.get_dof_indices(side_dof_indices);
1378 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1380 fill_n(local_rhs, ndofs, 0);
1381 local_flux_balance_vector.assign(ndofs, 0);
1382 local_flux_balance_rhs = 0;
1386 data_.bc_dirichlet_value[sbi].value_list(fe_values_side.
point_list(), ele_acc, bc_values);
1388 double side_flux = 0;
1389 for (
unsigned int k=0; k<qsize; k++)
1391 double transport_flux = side_flux/side.
measure();
1395 for (
unsigned int k=0; k<qsize; k++)
1397 double bc_term = -transport_flux*bc_values[k]*fe_values_side.
JxW(k);
1398 for (
unsigned int i=0; i<ndofs; i++)
1399 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k);
1401 for (
unsigned int i=0; i<ndofs; i++)
1402 local_flux_balance_rhs -= local_rhs[i];
1406 for (
unsigned int k=0; k<qsize; k++)
1408 double bc_term =
gamma[sbi][side.
cond_idx()]*bc_values[k]*fe_values_side.
JxW(k);
1410 for (
unsigned int i=0; i<ndofs; i++)
1411 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k)
1412 + arma::dot(bc_grad,fe_values_side.
shape_grad(i,k));
1414 for (
unsigned int k=0; k<qsize; k++)
1416 for (
unsigned int i=0; i<ndofs; i++)
1423 if (Model::time_->tlevel() > 0)
1424 for (
unsigned int i=0; i<ndofs; i++)
1425 local_flux_balance_rhs -= local_rhs[i];
1429 Model::get_flux_bc_data(sbi, fe_values_side.
point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1430 for (
unsigned int k=0; k<qsize; k++)
1432 double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.
JxW(k);
1433 for (
unsigned int i=0; i<ndofs; i++)
1434 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k);
1437 for (
unsigned int i=0; i<ndofs; i++)
1439 for (
unsigned int k=0; k<qsize; k++)
1440 local_flux_balance_vector[i] += csection[k]*bc_sigma[k]*fe_values_side.
JxW(k)*fe_values_side.
shape_value(i,k);
1441 local_flux_balance_rhs -= local_rhs[i];
1446 Model::get_flux_bc_data(sbi, fe_values_side.
point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1447 for (
unsigned int k=0; k<qsize; k++)
1449 double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.
JxW(k);
1450 for (
unsigned int i=0; i<ndofs; i++)
1451 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k);
1454 for (
unsigned int i=0; i<ndofs; i++)
1456 for (
unsigned int k=0; k<qsize; k++)
1457 local_flux_balance_vector[i] += csection[k]*(arma::dot(
ad_coef[sbi][k], fe_values_side.
normal_vector(k)) + bc_sigma[k])*fe_values_side.
JxW(k)*fe_values_side.
shape_value(i,k);
1458 local_flux_balance_rhs -= local_rhs[i];
1463 for (
unsigned int k=0; k<qsize; k++)
1465 for (
unsigned int i=0; i<ndofs; i++)
1471 Model::balance_->add_flux_matrix_values(Model::subst_idx[sbi], side, side_dof_indices, local_flux_balance_vector);
1472 Model::balance_->add_flux_vec_value(Model::subst_idx[sbi], side, local_flux_balance_rhs);
1480 template<
class Model>
1481 template<
unsigned int dim>
1486 OLD_ASSERT(cell->
dim() == dim,
"Element dimension mismatch!");
1490 for (
unsigned int k=0; k<fv.
n_points(); k++)
1492 velocity[k].zeros();
1493 for (
unsigned int sid=0; sid<cell->
n_sides(); sid++)
1494 for (
unsigned int c=0; c<3; ++c)
1501 template<
class Model>
1510 double delta = 0, h = 0;
1513 if (side.
dim() == 0)
1519 for (
unsigned int i=0; i<side.
n_nodes(); i++)
1520 for (
unsigned int j=i+1; j<side.
n_nodes(); j++)
1521 h = max(h, side.
node(i)->distance( *side.
node(j).node() ));
1525 for (
int k=0; k<K_size; k++)
1526 delta += dot(K[k]*normal_vector,normal_vector);
1536 template<
class Model>
1540 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1541 ls[sbi]->start_allocation();
1542 prepare_initial_condition<1>();
1543 prepare_initial_condition<2>();
1544 prepare_initial_condition<3>();
1546 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1547 ls[sbi]->start_add_assembly();
1548 prepare_initial_condition<1>();
1549 prepare_initial_condition<2>();
1550 prepare_initial_condition<3>();
1552 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1560 template<
class Model>
1561 template<
unsigned int dim>
1566 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
1568 double matrix[ndofs*ndofs],
rhs[ndofs];
1571 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1572 init_values[sbi].resize(qsize);
1574 for (
auto cell :
feo->
dh()->own_range() )
1576 if (cell.dim() != dim)
continue;
1579 cell.get_dof_indices(dof_indices);
1582 Model::compute_init_cond(fe_values.
point_list(), elem, init_values);
1584 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1586 for (
unsigned int i=0; i<ndofs; i++)
1589 for (
unsigned int j=0; j<ndofs; j++)
1590 matrix[i*ndofs+j] = 0;
1593 for (
unsigned int k=0; k<qsize; k++)
1595 double rhs_term = init_values[sbi][k]*fe_values.
JxW(k);
1597 for (
unsigned int i=0; i<ndofs; i++)
1599 for (
unsigned int j=0; j<ndofs; j++)
1605 ls[sbi]->
set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), matrix, rhs);
1611 template<
class Model>
1614 el_4_loc = Model::mesh_->get_el_4_loc();
1615 el_ds = Model::mesh_->get_el_ds();
1619 template<
class Model>
1622 if (solution_changed)
1624 unsigned int i_cell=0;
1625 for (
auto cell :
feo->
dh()->own_range() )
1628 unsigned int n_dofs;
1632 n_dofs =
feo->
fe<1>()->n_dofs();
1635 n_dofs =
feo->
fe<2>()->n_dofs();
1638 n_dofs =
feo->
fe<3>()->n_dofs();
1643 cell.get_dof_indices(dof_indices);
1645 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1647 double old_average = 0;
1648 for (
unsigned int j=0; j<n_dofs; ++j)
1649 old_average +=
ls[sbi]->get_solution_array()[dof_indices[j]-
feo->
dh()->distr()->begin()];
1650 old_average /= n_dofs;
1652 for (
unsigned int j=0; j<n_dofs; ++j)
1653 ls[sbi]->get_solution_array()[dof_indices[j]-
feo->
dh()->distr()->begin()] +=
solution_elem_[sbi][i_cell] - old_average;
1659 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1663 template<
class Model>
1666 return Model::mesh_->get_row_4_el();
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
const Edge * edge() const
Input::Record input_rec
Record with input specification.
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
static auto subdomain(Mesh &mesh) -> IndexField
void assemble_fluxes_element_element()
Assembles the fluxes between elements of the same dimension.
void set_sources()
Assembles the right hand side due to volume sources.
FiniteElement< dim > * fe()
void set_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
void get_par_info(LongIdx *&el_4_loc, Distribution *&el_ds)
unsigned int n_nodes() const
Transformed quadrature weight for cell sides.
vector< double > mm_coef
Mass matrix coefficients.
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.
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
unsigned int side_idx() const
void calculate_concentration_matrix()
Transport with dispersion implemented using discontinuous Galerkin method.
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
double distance(const Node &n2) const
void set_from_input(const Input::Record in_rec) override
virtual void start_add_assembly()
void assemble_mass_matrix()
Assembles the mass matrix.
void output(TimeStep step)
virtual PetscErrorCode mat_zero_entries()
vector< VectorMPI > output_vec
Array for storing the output solution data.
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
void update_solution() override
Computes the solution in one time instant.
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...
Fields computed from the mesh data.
virtual void start_allocation()
void set_initial_condition()
Sets the initial condition.
Cell accessor allow iterate over DOF handler cells.
Class FEValues calculates finite element data on the actual cells such as shape function values...
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
virtual void finish_assembly()=0
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
SideIter side(const unsigned int loc_index)
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
#define AVERAGE(i, k, side_id)
std::shared_ptr< DOFHandlerMultiDim > dh()
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
vector< vector< vector< arma::mat33 > > > dif_coef_edg
Diffusion coefficients on edges.
Discontinuous Galerkin method for equation of transport with dispersion.
const Vec & get_solution(unsigned int sbi)
bool is_own() const
Return true if accessor represents own element (false for ghost element)
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
ElementAccessor< 3 > element() const
vector< double > ret_sources_prev
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
static const Input::Type::Selection & get_dg_variant_selection_input_type()
Input type for the DG variant selection.
unsigned int n_points()
Returns the number of quadrature points.
EquationOutput output_fields
std::vector< Mat > mass_matrix
The mass matrix.
void calculate_cumulative_balance()
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
static constexpr bool value
Symmetric Gauss-Legendre quadrature formulae on simplices.
#define WAVERAGE(i, k, side_id)
TransportDG(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
vector< vector< double > > ret_coef
Retardation coefficient due to sorption.
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) override
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
void initialize() override
vector< vector< vector< arma::vec3 > > > ad_coef_edg
Advection coefficients on edges.
Transformed quadrature points.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
void output_data()
Postprocesses the solution and writes to output file.
virtual PetscErrorCode set_rhs(Vec &rhs)
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
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.
RangeConvert< DHEdgeSide, DHCellSide > edge_sides() const
Returns range of all sides looped over common Edge.
FieldCommon & input_default(const string &input_default)
const Vec & get_solution()
double JxW(const unsigned int point_no) override
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
Raviart-Thomas element of order 0.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Field< 3, FieldValue< 3 >::Scalar > subdomain
unsigned int n_sides() const
static auto region_id(Mesh &mesh) -> IndexField
#define START_TIMER(tag)
Starts a timer with specified tag.
std::vector< Vec > rhs
Vector of right hand side.
Field< 3, FieldValue< 3 >::Scalar > region_id
double elem_anisotropy(ElementAccessor< 3 > e)
unsigned int dg_order
Polynomial order of finite elements.
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Shape function gradients.
FLOW123D_FORCE_LINK_IN_CHILD(concentrationTransportModel)
NodeAccessor< 3 > node(unsigned int i) const
double measure() const
Calculate metrics of the side.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no) override
Returns the normal vector to a side at given quadrature point.
virtual PetscErrorCode rhs_zero_entries()
void set_solution(Vec sol_vec)
Discontinuous Galerkin method for equation of transport with dispersion.
double ** solution_elem_
Element averages of solution (the array is passed to reactions in operator splitting).
std::vector< Mat > stiffness_matrix
The stiffness matrix.
FieldCommon & description(const string &description)
void set_values(int nrow, int *rows, int ncol, int *cols, PetscScalar *mat_vals, PetscScalar *rhs_vals)
Set values in the system matrix and values in the right-hand side vector on corresponding rows...
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.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
Discontinuous Galerkin method for equation of transport with dispersion.
double shape_value(const unsigned int function_no, const unsigned int point_no) override
Return the value of the function_no-th shape function at the point_no-th quadrature point...
void update_after_reactions(bool solution_changed)
MappingP1< dim, 3 > * mapping()
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the value of the function_no-th shape function at the point_no-th quadrature point...
vector< double > ret_sources
Temporary values of increments due to retardation (e.g. sorption)
bool set_time(const TimeStep &time, LimitSide limit_side)
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining 'warning' record of log.
FieldCommon & name(const string &name)
virtual SolveInfo solve()=0
#define END_TIMER(tag)
Ends a timer with specified tag.
Discontinuous Galerkin method for equation of transport with dispersion.
void calculate_velocity(const ElementAccessor< 3 > &cell, std::vector< arma::vec3 > &velocity, FEValuesBase< dim, 3 > &fv)
Calculates the velocity field on a given dim dimensional cell.
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
void set_DG_parameters_boundary(Side side, const int K_size, const std::vector< arma::mat33 > &K, const double flux, const arma::vec3 &normal_vector, const double alpha, double &gamma)
Calculates the dispersivity (diffusivity) tensor from the velocity field.
static const Input::Type::Record & get_input_type()
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
std::vector< Vec > mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time)...
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Base class for FEValues and FESideValues.
static UnitSI & dimensionless()
Returns dimensionless unit.
static bool print_message_table(ostream &stream, std::string equation_name)
~TransportDG()
Destructor.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
LinSys ** ls
Linear algebra system for the transport equation.
unsigned int get_dof_indices(std::vector< int > &indices) const
Fill vector of the global indices of dofs associated to the cell.
FiniteElement< dim > * fe_rt()
SideIter side(const unsigned int i) const
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
#define JUMP(i, k, side_id)
Definitions of Raviart-Thomas finite elements.
Side accessor allows to iterate over sides of DOF handler cell.
ElementAccessor< 3 > element_accessor()
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
unsigned int n_nodes() const
Transformed quadrature weights.
Calculates finite element data on a side.
const Node * node(unsigned int ni) const