49 return 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")
75 std::string equation_name = std::string(Model::ModelEqData::name()) +
"_DG";
76 return Model::get_input_type(
"DG",
"Discontinuous Galerkin (DG) solver")
78 "Solver for the linear system.")
81 .make_field_descriptor_type(equation_name)),
83 "Input fields of the equation.")
85 "Variant of the interior penalty discontinuous Galerkin method.")
87 "Polynomial order for the finite element in DG method (order 0 is suitable if there is no diffusion/dispersion).")
89 EqData().output_fields.make_output_type(equation_name,
""),
90 IT::Default(
"{ \"fields\": [ " + Model::ModelEqData::default_output_field() +
"] }"),
91 "Specification of output fields and output times.")
97 Input::register_class< TransportDG<Model>,
Mesh &,
const Input::Record>(std::string(Model::ModelEqData::name()) +
"_DG") +
105 unsigned int q_order;
139 xprintf(
PrgErr,
"Unsupported polynomial order %d for finite elements in TransportDG ", fe_order);
159 dh_->distribute_dofs(*fe1_, *fe2_, *fe3_);
205 template<
class Model>
209 .
name(
"fracture_sigma")
211 "Coefficient of diffusive transfer through fractures (for each substance).")
219 "Penalty parameter influencing the discontinuity of the solution (for each substance). " 220 "Its default value 1 is sufficient in most cases. Higher value diminishes the inter-element jumps.")
233 .
description(
"Subdomain ids of the domain decomposition.");
241 template<
class Model>
243 : Model(init_mesh, in_rec),
257 data_.set_mesh(init_mesh);
266 Model::init_from_input(in_rec);
275 template<
class Model>
278 data_.set_components(Model::substances_.names());
282 gamma.resize(Model::n_substances());
283 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
284 gamma[sbi].resize(Model::mesh_->boundary_.size());
287 int qsize = max(
feo->
q<0>()->size(), max(
feo->
q<1>()->size(), max(
feo->
q<2>()->size(),
feo->
q<3>()->size())));
288 int max_edg_sides = max(Model::mesh_->max_edge_sides(1), max(Model::mesh_->max_edge_sides(2), Model::mesh_->max_edge_sides(3)));
290 ret_coef.resize(Model::n_substances());
293 ad_coef.resize(Model::n_substances());
294 dif_coef.resize(Model::n_substances());
295 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
303 for (
int sd=0; sd<max_edg_sides; sd++)
307 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
319 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
323 VecCreateSeq(PETSC_COMM_SELF, output_vector_size, &
output_vec[sbi]);
325 data_.output_field.set_components(Model::substances_.names());
326 data_.output_field.set_mesh(*Model::mesh_);
329 data_.output_field.setup_components();
330 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
335 data_.output_field[sbi].set_field(Model::mesh_->region_db().get_region_set(
"ALL"), output_field_ptr, 0);
342 std::string petsc_default_opts;
344 petsc_default_opts =
"-ksp_type bcgs -pc_type ilu -pc_factor_levels 2 -ksp_diagonal_scale_fix -pc_factor_fill 6.0";
346 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";
349 ls =
new LinSys*[Model::n_substances()];
354 mass_matrix.resize(Model::n_substances(),
nullptr);
355 rhs.resize(Model::n_substances(),
nullptr);
356 mass_vec.resize(Model::n_substances(),
nullptr);
357 ret_vec.resize(Model::n_substances(),
nullptr);
359 for (
unsigned int sbi = 0; sbi < Model::n_substances(); sbi++) {
366 solution_elem_[sbi] =
new double[Model::mesh_->get_el_ds()->lsize()];
374 max(
feo->
fe<1>()->n_dofs(), max(
feo->
fe<2>()->n_dofs(),
feo->
fe<3>()->n_dofs())));
379 template<
class Model>
384 if (
gamma.size() > 0) {
387 for (
auto &vec :
output_vec) VecDestroy(&vec);
389 for (
unsigned int i=0; i<Model::n_substances(); i++)
421 template<
class Model>
424 VecScatter output_scatter;
425 VecScatterCreateToZero(
ls[0]->
get_solution(), &output_scatter, PETSC_NULL);
426 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
432 VecScatterDestroy(&(output_scatter));
438 template<
class Model>
442 data_.mark_input_times( *(Model::time_) );
444 std::stringstream ss;
452 for (
unsigned int sbi = 0; sbi < Model::n_substances(); sbi++)
459 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
461 Model::balance_->calculate_instant(Model::subst_idx[sbi],
ls[sbi]->
get_solution());
471 template<
class Model>
475 for (
unsigned int i=0; i<Model::n_substances(); i++)
491 for (
unsigned int i=0; i<Model::n_substances(); i++)
502 template<
class Model>
507 Model::time_->next_time();
508 Model::time_->view(
"TDG");
517 for (
unsigned int i=0; i<Model::n_substances(); i++)
524 for (
unsigned int i=0; i<Model::n_substances(); i++)
534 MatConvert(*(
ls_dt[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &
mass_matrix[i]);
537 MatCopy(*(
ls_dt[i]->get_matrix() ),
mass_matrix[i], DIFFERENT_NONZERO_PATTERN);
544 || Model::flux_changed)
548 for (
unsigned int i=0; i<Model::n_substances(); i++)
554 for (
unsigned int i=0; i<Model::n_substances(); i++)
559 MatConvert(*(
ls[i]->get_matrix() ), MATSAME, MAT_INITIAL_MATRIX, &
stiffness_matrix[i]);
568 || Model::flux_changed)
570 for (
unsigned int i=0; i<Model::n_substances(); i++)
577 for (
unsigned int i=0; i<Model::n_substances(); i++)
581 if (
rhs[i] ==
nullptr) VecDuplicate(*(
ls[i]->get_rhs() ), &
rhs[i]);
582 VecCopy(*(
ls[i]->get_rhs() ),
rhs[i]);
586 Model::flux_changed =
false;
607 for (
unsigned int i=0; i<Model::n_substances(); i++)
610 MatAXPY(m, 1./Model::time_->dt(),
mass_matrix[i], SUBSET_NONZERO_PATTERN);
613 VecDuplicate(
rhs[i], &w);
614 VecWAXPY(w, 1./Model::time_->dt(),
mass_vec[i],
rhs[i]);
633 template<
class Model>
637 for (
unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
645 n_dofs =
feo->
fe<1>()->n_dofs();
648 n_dofs =
feo->
fe<2>()->n_dofs();
651 n_dofs =
feo->
fe<3>()->n_dofs();
655 unsigned int dof_indices[n_dofs];
658 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
662 for (
unsigned int j=0; j<n_dofs; ++j)
673 template<
class Model>
687 Model::output_data();
690 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
691 Model::balance_->calculate_instant(Model::subst_idx[sbi],
ls[sbi]->
get_solution());
692 Model::balance_->output();
699 template<
class Model>
702 if (Model::balance_->cumulative())
704 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
706 Model::balance_->calculate_cumulative(Model::subst_idx[sbi],
ls[sbi]->
get_solution());
718 template<
class Model>
722 Model::balance_->start_mass_assembly(Model::subst_idx);
723 assemble_mass_matrix<1>();
724 assemble_mass_matrix<2>();
725 assemble_mass_matrix<3>();
726 Model::balance_->finish_mass_assembly(Model::subst_idx);
731 template<
class Model>
template<
unsigned int dim>
735 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
737 PetscScalar local_mass_matrix[ndofs*ndofs], local_retardation_balance_vector[ndofs];
741 for (
unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
744 if (cell->dim() != dim)
continue;
750 Model::compute_mass_matrix_coefficient(fe_values.
point_list(), ele_acc,
mm_coef);
753 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
756 for (
unsigned int i=0; i<ndofs; i++)
758 for (
unsigned int j=0; j<ndofs; j++)
760 local_mass_matrix[i*ndofs+j] = 0;
761 for (
unsigned int k=0; k<qsize; k++)
766 for (
unsigned int i=0; i<ndofs; i++)
768 local_mass_balance_vector[i] = 0;
769 local_retardation_balance_vector[i] = 0;
770 for (
unsigned int k=0; k<qsize; k++)
777 Model::balance_->add_mass_matrix_values(Model::subst_idx[sbi], ele_acc.
region().
bulk_idx(), dof_indices, local_mass_balance_vector);
778 ls_dt[sbi]->
mat_set_values(ndofs, &(dof_indices[0]), ndofs, &(dof_indices[0]), local_mass_matrix);
779 VecSetValues(
ret_vec[sbi], ndofs, &(dof_indices[0]), local_retardation_balance_vector, ADD_VALUES);
787 template<
class Model>
792 assemble_volume_integrals<1>();
793 assemble_volume_integrals<2>();
794 assemble_volume_integrals<3>();
798 assemble_fluxes_boundary<1>();
799 assemble_fluxes_boundary<2>();
800 assemble_fluxes_boundary<3>();
804 assemble_fluxes_element_element<1>();
805 assemble_fluxes_element_element<2>();
806 assemble_fluxes_element_element<3>();
810 assemble_fluxes_element_side<1>();
811 assemble_fluxes_element_side<2>();
812 assemble_fluxes_element_side<3>();
819 template<
class Model>
820 template<
unsigned int dim>
827 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
828 unsigned int dof_indices[ndofs];
831 PetscScalar local_matrix[ndofs*ndofs];
834 for (
unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
837 if (cell->dim() != dim)
continue;
846 Model::compute_sources_sigma(fe_values.
point_list(), ele_acc, sources_sigma);
849 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
851 for (
unsigned int i=0; i<ndofs; i++)
852 for (
unsigned int j=0; j<ndofs; j++)
853 local_matrix[i*ndofs+j] = 0;
855 for (
unsigned int k=0; k<qsize; k++)
857 for (
unsigned int i=0; i<ndofs; i++)
860 double ad_dot_grad_i = arma::dot(ad_coef[sbi][k], fe_values.
shape_grad(i,k));
862 for (
unsigned int j=0; j<ndofs; j++)
863 local_matrix[i*ndofs+j] += (arma::dot(Kt_grad_i, fe_values.
shape_grad(j,k))
868 ls[sbi]->
mat_set_values(ndofs, (
int *)dof_indices, ndofs, (
int *)dof_indices, local_matrix);
874 template<
class Model>
878 Model::balance_->start_source_assembly(Model::subst_idx);
882 Model::balance_->finish_source_assembly(Model::subst_idx);
886 template<
class Model>
887 template<
unsigned int dim>
892 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
897 PetscScalar local_rhs[ndofs];
902 for (
unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
905 if (cell->dim() != dim)
continue;
910 Model::compute_source_coefficients(fe_values.
point_list(), cell->element_accessor(), sources_conc, sources_density, sources_sigma);
913 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
915 fill_n(local_rhs, ndofs, 0);
916 local_source_balance_vector.assign(ndofs, 0);
917 local_source_balance_rhs.assign(ndofs, 0);
920 for (
unsigned int k=0; k<qsize; k++)
922 source = (sources_density[sbi][k] + sources_conc[sbi][k]*sources_sigma[sbi][k])*fe_values.
JxW(k);
924 for (
unsigned int i=0; i<ndofs; i++)
929 for (
unsigned int i=0; i<ndofs; i++)
931 for (
unsigned int k=0; k<qsize; k++)
932 local_source_balance_vector[i] -= sources_sigma[sbi][k]*fe_values.
shape_value(i,k)*fe_values.
JxW(k);
934 local_source_balance_rhs[i] += local_rhs[i];
936 Model::balance_->add_source_matrix_values(Model::subst_idx[sbi], cell->region().bulk_idx(), dof_indices, local_source_balance_vector);
937 Model::balance_->add_source_vec_values(Model::subst_idx[sbi], cell->region().bulk_idx(), dof_indices, local_source_balance_rhs);
944 template<
class Model>
945 template<
unsigned int dim>
951 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim-1>()->size(),
954 PetscScalar local_matrix[ndofs*ndofs];
957 double gamma_l, omega[2], transport_flux;
959 for (
unsigned int sid=0; sid<n_max_sides; sid++)
961 side_dof_indices.push_back(
new unsigned int[ndofs]);
972 for (
int sid=0; sid<edg->
n_sides; sid++)
977 fe_values[sid]->reinit(cell, edg->
side(sid)->
el_idx());
980 Model::compute_advection_diffusion_coefficients(fe_values[sid]->point_list(), side_velocity[sid], ele_acc,
ad_coef_edg[sid],
dif_coef_edg[sid]);
981 dg_penalty[sid].resize(Model::n_substances());
982 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
983 dg_penalty[sid][sbi] =
data_.
dg_penalty[sbi].value(cell->centre(), ele_acc);
987 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
990 for (
int sid=0; sid<edg->
n_sides; sid++)
993 for (
unsigned int k=0; k<qsize; k++)
994 fluxes[sid] += arma::dot(
ad_coef_edg[sid][sbi][k], fe_values[sid]->normal_vector(k))*fe_values[sid]->JxW(k);
998 for (
int s1=0; s1<edg->
n_sides; s1++)
1000 for (
int s2=s1+1; s2<edg->
n_sides; s2++)
1006 arma::vec3 nv = fe_values[s1]->normal_vector(0);
1009 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);
1015 #define AVERAGE(i,k,side_id) (fe_values[sd[side_id]]->shape_value(i,k)*0.5) 1016 #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]) 1017 #define JUMP(i,k,side_id) ((side_id==0?1:-1)*fe_values[sd[side_id]]->shape_value(i,k)) 1020 for (
int n=0; n<2; n++)
1024 for (
int m=0; m<2; m++)
1026 for (
unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
1027 for (
unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
1028 local_matrix[i*fe_values[sd[m]]->n_dofs()+j] = 0;
1030 for (
unsigned int k=0; k<qsize; k++)
1032 double flux_times_JxW = transport_flux*fe_values[0]->JxW(k);
1033 double gamma_times_JxW = gamma_l*fe_values[0]->JxW(k);
1035 for (
unsigned int i=0; i<fe_values[sd[n]]->n_dofs(); i++)
1037 double flux_JxW_jump_i = flux_times_JxW*
JUMP(i,k,n);
1038 double gamma_JxW_jump_i = gamma_times_JxW*
JUMP(i,k,n);
1039 double JxW_jump_i = fe_values[0]->JxW(k)*
JUMP(i,k,n);
1042 for (
unsigned int j=0; j<fe_values[sd[m]]->n_dofs(); j++)
1044 int index = i*fe_values[sd[m]]->n_dofs()+j;
1047 local_matrix[index] += flux_JxW_jump_i*
AVERAGE(j,k,m);
1050 local_matrix[index] += gamma_JxW_jump_i*
JUMP(j,k,m);
1053 local_matrix[index] -=
WAVERAGE(j,k,m)*JxW_jump_i;
1054 local_matrix[index] -=
JUMP(j,k,m)*JxW_var_wavg_i;
1058 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);
1069 for (
unsigned int i=0; i<n_max_sides; i++)
1071 delete fe_values[i];
1072 delete[] side_dof_indices[i];
1077 template<
class Model>
1078 template<
unsigned int dim>
1085 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim-1>()->size();
1086 unsigned int side_dof_indices[ndofs];
1087 PetscScalar local_matrix[ndofs*ndofs];
1097 if (edg->
n_sides > 1)
continue;
1099 if (edg->
side(0)->
dim() != dim-1)
continue;
1101 if (edg->
side(0)->
cond() == NULL)
continue;
1111 Model::compute_advection_diffusion_coefficients(fe_values_side.
point_list(), side_velocity, ele_acc,
ad_coef,
dif_coef);
1114 data_.cross_section.value_list(fe_values_side.
point_list(), ele_acc, csection);
1116 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1118 for (
unsigned int i=0; i<ndofs; i++)
1119 for (
unsigned int j=0; j<ndofs; j++)
1120 local_matrix[i*ndofs+j] = 0;
1124 double side_flux = 0;
1125 for (
unsigned int k=0; k<qsize; k++)
1126 side_flux += arma::dot(ad_coef[sbi][k], fe_values_side.
normal_vector(k))*fe_values_side.
JxW(k);
1127 double transport_flux = side_flux/side->
measure();
1134 transport_flux += gamma_l;
1138 for (
unsigned int k=0; k<qsize; k++)
1140 double flux_times_JxW;
1144 flux_times_JxW = csection[k]*robin_sigma[k]*fe_values_side.
JxW(k);
1149 flux_times_JxW = (transport_flux + csection[k]*robin_sigma[k])*fe_values_side.
JxW(k);
1154 flux_times_JxW = transport_flux*fe_values_side.
JxW(k);
1156 for (
unsigned int i=0; i<ndofs; i++)
1158 for (
unsigned int j=0; j<ndofs; j++)
1161 local_matrix[i*ndofs+j] += flux_times_JxW*fe_values_side.
shape_value(i,k)*fe_values_side.
shape_value(j,k);
1167 )*fe_values_side.
JxW(k);
1172 ls[sbi]->
mat_set_values(ndofs, (
int *)side_dof_indices, ndofs, (
int *)side_dof_indices, local_matrix);
1178 template<
class Model>
1179 template<
unsigned int dim>
1183 if (dim == 1)
return;
1194 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs();
1195 const unsigned int qsize =
feo->
q<dim-1>()->size();
1196 unsigned int side_dof_indices[2*ndofs], n_dofs[2];
1200 PetscScalar local_matrix[4*ndofs*ndofs];
1201 double comm_flux[2][2];
1205 fv_sb[0] = &fe_values_vb;
1206 fv_sb[1] = &fe_values_side;
1217 fe_values_vb.
reinit(cell_sub);
1218 n_dofs[0] = fv_sb[0]->n_dofs();
1223 n_dofs[1] = fv_sb[1]->n_dofs();
1227 element_id[0] = cell_sub.
index();
1228 element_id[1] = cell.
index();
1236 data_.cross_section.value_list(fe_values_vb.
point_list(), cell_sub->element_accessor(), csection_lower);
1237 data_.cross_section.value_list(fe_values_vb.
point_list(), cell->element_accessor(), csection_higher);
1239 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1241 for (
unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
1242 for (
unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
1243 local_matrix[i*(n_dofs[0]+n_dofs[1])+j] = 0;
1248 for (
unsigned int k=0; k<qsize; k++)
1259 2*csection_higher[k]*csection_higher[k]/(csection_lower[k]*csection_lower[k]);
1263 comm_flux[0][0] = (sigma-min(0.,transport_flux))*fv_sb[0]->JxW(k);
1264 comm_flux[0][1] = -(sigma-min(0.,transport_flux))*fv_sb[0]->JxW(k);
1265 comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1266 comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb[0]->JxW(k);
1268 for (
int n=0; n<2; n++)
1272 for (
unsigned int i=0; i<n_dofs[n]; i++)
1273 for (
int m=0; m<2; m++)
1274 for (
unsigned int j=0; j<n_dofs[m]; j++)
1275 local_matrix[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
1276 comm_flux[m][n]*fv_sb[m]->shape_value(j,k)*fv_sb[n]->shape_value(i,k);
1279 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);
1290 template<
class Model>
1294 Model::balance_->start_flux_assembly(Model::subst_idx);
1295 set_boundary_conditions<1>();
1296 set_boundary_conditions<2>();
1297 set_boundary_conditions<3>();
1298 Model::balance_->finish_flux_assembly(Model::subst_idx);
1303 template<
class Model>
1304 template<
unsigned int dim>
1311 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim-1>()->size();
1313 unsigned int loc_b=0;
1314 double local_rhs[ndofs];
1316 PetscScalar local_flux_balance_rhs;
1320 bc_ref_values(qsize);
1324 for (
unsigned int loc_el = 0; loc_el < Model::mesh_->get_el_ds()->lsize(); loc_el++)
1327 if (elm->boundary_idx_ ==
nullptr)
continue;
1331 Edge *edg = elm->side(si)->edge();
1332 if (edg->
n_sides > 1)
continue;
1334 if (edg->
side(0)->
cond() == NULL)
continue;
1336 if (edg->
side(0)->
dim() != dim-1)
1338 if (edg->
side(0)->
cond() !=
nullptr) ++loc_b;
1347 Model::get_bc_type(ele_acc, bc_type);
1354 data_.cross_section.value_list(fe_values_side.
point_list(), side->
element()->element_accessor(), csection);
1358 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1360 fill_n(local_rhs, ndofs, 0);
1361 local_flux_balance_vector.assign(ndofs, 0);
1362 local_flux_balance_rhs = 0;
1366 data_.bc_dirichlet_value[sbi].value_list(fe_values_side.
point_list(), ele_acc, bc_values);
1368 double side_flux = 0;
1369 for (
unsigned int k=0; k<qsize; k++)
1371 double transport_flux = side_flux/side->
measure();
1375 for (
unsigned int k=0; k<qsize; k++)
1377 double bc_term = -transport_flux*bc_values[k]*fe_values_side.
JxW(k);
1378 for (
unsigned int i=0; i<ndofs; i++)
1379 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k);
1381 for (
unsigned int i=0; i<ndofs; i++)
1382 local_flux_balance_rhs -= local_rhs[i];
1386 for (
unsigned int k=0; k<qsize; k++)
1388 double bc_term =
gamma[sbi][side->
cond_idx()]*bc_values[k]*fe_values_side.
JxW(k);
1390 for (
unsigned int i=0; i<ndofs; i++)
1391 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k)
1392 + arma::dot(bc_grad,fe_values_side.
shape_grad(i,k));
1394 for (
unsigned int k=0; k<qsize; k++)
1396 for (
unsigned int i=0; i<ndofs; i++)
1403 if (Model::time_->tlevel() > 0)
1404 for (
unsigned int i=0; i<ndofs; i++)
1405 local_flux_balance_rhs -= local_rhs[i];
1409 Model::get_flux_bc_data(sbi, fe_values_side.
point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1410 for (
unsigned int k=0; k<qsize; k++)
1412 double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.
JxW(k);
1413 for (
unsigned int i=0; i<ndofs; i++)
1414 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k);
1417 for (
unsigned int i=0; i<ndofs; i++)
1419 for (
unsigned int k=0; k<qsize; k++)
1420 local_flux_balance_vector[i] += csection[k]*bc_sigma[k]*fe_values_side.
JxW(k)*fe_values_side.
shape_value(i,k);
1421 local_flux_balance_rhs -= local_rhs[i];
1426 Model::get_flux_bc_data(sbi, fe_values_side.
point_list(), ele_acc, bc_fluxes, bc_sigma, bc_ref_values);
1427 for (
unsigned int k=0; k<qsize; k++)
1429 double bc_term = csection[k]*(bc_sigma[k]*bc_ref_values[k]+bc_fluxes[k])*fe_values_side.
JxW(k);
1430 for (
unsigned int i=0; i<ndofs; i++)
1431 local_rhs[i] += bc_term*fe_values_side.
shape_value(i,k);
1434 for (
unsigned int i=0; i<ndofs; i++)
1436 for (
unsigned int k=0; k<qsize; k++)
1437 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);
1438 local_flux_balance_rhs -= local_rhs[i];
1443 for (
unsigned int k=0; k<qsize; k++)
1445 for (
unsigned int i=0; i<ndofs; i++)
1451 Model::balance_->add_flux_matrix_values(Model::subst_idx[sbi], loc_b, side_dof_indices, local_flux_balance_vector);
1452 Model::balance_->add_flux_vec_value(Model::subst_idx[sbi], loc_b, local_flux_balance_rhs);
1461 template<
class Model>
1462 template<
unsigned int dim>
1465 OLD_ASSERT(cell->dim() == dim,
"Element dimension mismatch!");
1469 for (
unsigned int k=0; k<fv.
n_points(); k++)
1471 velocity[k].zeros();
1472 for (
unsigned int sid=0; sid<cell->n_sides(); sid++)
1473 velocity[k] += fv.
shape_vector(sid,k) * Model::mh_dh->side_flux( *(cell->side(sid)) );
1481 double h_max = 0, h_min = numeric_limits<double>::infinity();
1482 for (
unsigned int i=0; i<e->n_nodes(); i++)
1483 for (
unsigned int j=i+1; j<e->n_nodes(); j++)
1485 h_max = max(h_max, e->node[i]->distance(*e->node[j]));
1486 h_min = min(h_min, e->node[i]->distance(*e->node[j]));
1493 template<
class Model>
1502 const double alpha1,
1503 const double alpha2,
1506 double &transport_flux)
1510 double local_alpha = 0;
1522 for (
unsigned int i=0; i<s->n_nodes(); i++)
1523 for (
unsigned int j=i+1; j<s->n_nodes(); j++)
1524 h = max(h, s->node(i)->distance(*s->node(j)));
1532 double pflux = 0, nflux = 0;
1533 for (
int i=0; i<edg.
n_sides; i++)
1542 if (fluxes[s2] > 0 && fluxes[s1] < 0 && s1 < s2)
1543 transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
1544 else if (fluxes[s2] < 0 && fluxes[s1] > 0 && s1 < s2)
1545 transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
1547 transport_flux = fluxes[s1];
1551 gamma = 0.5*fabs(transport_flux);
1555 local_alpha = max(alpha1, alpha2);
1563 for (
int k=0; k<K_size; k++)
1564 delta[0] += dot(K1[k]*normal_vector,normal_vector);
1567 gamma += local_alpha/h*aniso1*aniso2*delta[0];
1573 for (
int k=0; k<K_size; k++)
1575 delta[0] += dot(K1[k]*normal_vector,normal_vector);
1576 delta[1] += dot(K2[k]*normal_vector,normal_vector);
1581 double delta_sum = delta[0] + delta[1];
1584 if (fabs(delta_sum) > 0)
1586 omega[0] = delta[1]/delta_sum;
1587 omega[1] = delta[0]/delta_sum;
1588 gamma += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
1591 for (
int i=0; i<2; i++) omega[i] = 0;
1600 template<
class Model>
1609 double delta = 0, h = 0;
1612 if (side->
dim() == 0)
1618 for (
unsigned int i=0; i<side->
n_nodes(); i++)
1619 for (
unsigned int j=i+1; j<side->
n_nodes(); j++)
1624 for (
int k=0; k<K_size; k++)
1625 delta += dot(K[k]*normal_vector,normal_vector);
1635 template<
class Model>
1639 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1640 ls[sbi]->start_allocation();
1641 prepare_initial_condition<1>();
1642 prepare_initial_condition<2>();
1643 prepare_initial_condition<3>();
1645 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1646 ls[sbi]->start_add_assembly();
1647 prepare_initial_condition<1>();
1648 prepare_initial_condition<2>();
1649 prepare_initial_condition<3>();
1651 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1659 template<
class Model>
1660 template<
unsigned int dim>
1665 const unsigned int ndofs =
feo->
fe<dim>()->n_dofs(), qsize =
feo->
q<dim>()->size();
1666 unsigned int dof_indices[ndofs];
1667 double matrix[ndofs*ndofs],
rhs[ndofs];
1670 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1671 init_values[sbi].resize(qsize);
1673 for (
unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
1676 if (elem->dim() != dim)
continue;
1682 Model::compute_init_cond(fe_values.
point_list(), ele_acc, init_values);
1684 for (
unsigned int sbi=0; sbi<Model::n_substances(); sbi++)
1686 for (
unsigned int i=0; i<ndofs; i++)
1689 for (
unsigned int j=0; j<ndofs; j++)
1690 matrix[i*ndofs+j] = 0;
1693 for (
unsigned int k=0; k<qsize; k++)
1695 double rhs_term = init_values[sbi][k]*fe_values.
JxW(k);
1697 for (
unsigned int i=0; i<ndofs; i++)
1699 for (
unsigned int j=0; j<ndofs; j++)
1705 ls[sbi]->
set_values(ndofs, (
int *)dof_indices, ndofs, (
int *)dof_indices, matrix, rhs);
1711 template<
class Model>
1714 el_4_loc = Model::mesh_->get_el_4_loc();
1715 el_ds = Model::mesh_->get_el_ds();
1719 template<
class Model>
1722 if (solution_changed)
1724 for (
unsigned int i_cell=0; i_cell<Model::mesh_->get_el_ds()->lsize(); i_cell++)
1728 unsigned int n_dofs;
1729 switch (elem->dim())
1732 n_dofs =
feo->
fe<1>()->n_dofs();
1735 n_dofs =
feo->
fe<2>()->n_dofs();
1738 n_dofs =
feo->
fe<3>()->n_dofs();
1742 unsigned int dof_indices[n_dofs];
1745 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1747 double old_average = 0;
1748 for (
unsigned int j=0; j<n_dofs; ++j)
1749 old_average +=
ls[sbi]->get_solution_array()[dof_indices[j]-
feo->
dh()->
distr()->
begin()];
1750 old_average /= n_dofs;
1752 for (
unsigned int j=0; j<n_dofs; ++j)
1758 for (
unsigned int sbi=0; sbi<Model::n_substances(); ++sbi)
1762 template<
class Model>
1765 return Model::mesh_->get_row_4_el();
Input::Record input_rec
Record with input specification.
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
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.
void set_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
Transformed quadrature weight for cell sides.
FiniteElement< dim, 3 > * fe()
vector< double > mm_coef
Mass matrix coefficients.
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...
void calculate_concentration_matrix()
Transport with dispersion implemented using discontinuous Galerkin method.
Field< 3, FieldValue< 3 >::Integer > region_id
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()
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.
void initialize(std::shared_ptr< OutputTime > stream, Input::Record in_rec, const TimeGovernor &tg)
Class FEValues calculates finite element data on the actual cells such as shape function values...
#define FOR_ELEMENT_SIDES(i, j)
virtual void finish_assembly()=0
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
unsigned int n_loc_nb() const
Returns number of local neighbours.
#define AVERAGE(i, k, side_id)
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
int el_index(int loc_el) const
Returns the global index of local element.
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)
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
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 > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
void initialize() override
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.
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.
FieldCommon & input_default(const string &input_default)
const Vec & get_solution()
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.
unsigned int begin(int proc) const
get starting local index
int nb_index(int loc_nb) const
Returns the global index of local neighbour.
static auto region_id(Mesh &mesh) -> IndexField
#define START_TIMER(tag)
Starts a timer with specified tag.
Provides the numbering of the finite element degrees of freedom on the computational mesh...
std::vector< Vec > rhs
Vector of right hand side.
FiniteElement< dim, 3 > * fe_rt()
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.
Affine mapping between reference and actual cell.
unsigned int np() const
get num of processors
Shape function gradients.
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...
FLOW123D_FORCE_LINK_IN_CHILD(concentrationTransportModel)
void set_solution(double *sol_array)
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...
double elem_anisotropy(const ElementFullIter &e)
DOFHandlerMultiDim * dh()
virtual PetscErrorCode rhs_zero_entries()
unsigned int n_loc_edges() const
Returns number of local edges.
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.
Discontinuous Galerkin method for equation of transport with dispersion.
void get_dof_indices(const CellIterator &cell, unsigned int indices[]) const override
Returns the global indices of dofs associated to the cell.
ElementFullIter element() const
void update_after_reactions(bool solution_changed)
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
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.
Abstract linear system class.
Definitions of particular quadrature rules on simplices.
#define WarningOut()
Macro defining 'warning' record of log.
FieldCommon & name(const string &name)
vector< Vec > output_vec
Array for storing the output solution data.
#define END_TIMER(tag)
Ends a timer with specified tag.
Discontinuous Galerkin method for equation of transport with dispersion.
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
static const Input::Type::Record & get_input_type()
unsigned int el_idx() const
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)
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
Field< 3, FieldValue< 3 >::Integer > subdomain
bool el_is_local(int index) const
int edge_index(int loc_edg) const
Returns the global index of local edge.
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.
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...
Distribution * distr() const
Mapping< dim, 3 > * mapping()
static UnitSI & dimensionless()
Returns dimensionless unit.
static bool print_message_table(ostream &stream, std::string equation_name)
~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)
Calculates the dispersivity (diffusivity) tensor from the velocity field.
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)
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.
void get_par_info(int *&el_4_loc, Distribution *&el_ds)
const unsigned int n_global_dofs() const
Getter for the number of all mesh dofs required by the given finite element.
ElementAccessor< 3 > element_accessor()
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)
unsigned int n_nodes() const
Transformed quadrature weights.
unsigned int lsize(int proc) const
get local size