34 #include "petscviewer.h"
35 #include "petscerror.h"
37 #include <boost/foreach.hpp>
82 namespace it = Input::Type;
86 .
add_value(NoMortar,
"None",
"Mortar space: P0 on elements of lower dimension.")
87 .
add_value(MortarP0,
"P0",
"Mortar space: P0 on elements of lower dimension.")
88 .
add_value(MortarP1,
"P1",
"Mortar space: P1 on intersections, using non-conforming pressures.");
103 "Number of Schur complements to perform when solving MH sytem.")
105 "Linear solver for MH problem.")
107 "Parameters of output form MH module.")
109 "Method for coupling Darcy flow between dimensions." );
113 =
it::Record(
"Steady_MH",
"Mixed-Hybrid solver for STEADY saturated Darcy flow.")
123 =
it::Record(
"Unsteady_MH",
"Mixed-Hybrid solver for unsteady saturated Darcy flow.")
126 "Time governor setting for the unsteady Darcy flow model.")
127 .
copy_keys(DarcyFlowMH_Steady::input_type);
131 =
it::Record(
"Unsteady_LMH",
"Lumped Mixed-Hybrid solver for unsteady saturated Darcy flow.")
134 "Time governor setting for the unsteady Darcy flow model.")
135 .
copy_keys(DarcyFlowMH_Steady::input_type);
149 ADD_FIELD(
cross_section,
"Complement dimension parameter (cross section for 1D, thickness for 2D).",
"1.0" );
151 ADD_FIELD(
sigma,
"Transition coefficient between dimensions.",
"1.0");
154 ADD_FIELD(
bc_type,
"Boundary condition type, possible values:",
"\"none\"" );
176 main_matrix_fields = this->subset({
"anisotropy",
"conductivity",
"cross_section",
"sigma",
"bc_type",
"bc_robin_sigma"});
177 rhs_fields = this->subset({
"water_source_density",
"bc_pressure",
"bc_flux"});
200 using namespace Input;
277 xprintf(
MsgLog,
"Linear solver ended with reason: %d \n", convergedReason );
278 ASSERT( convergedReason >= 0,
"Linear solver failed to converge. Convergence reason %d \n", convergedReason );
340 schur0 -> get_whole_solution( sol_disordered );
344 for (
int i = 0; i < this->
size; i++ ) {
351 ASSERT(vec != NULL,
"Requested solution is not allocated!\n");
365 ASSERT(vec != NULL,
"Requested solution is not allocated!\n");
411 int el_row, side_row, edge_row;
414 int side_rows[4], edge_rows[4];
417 double minus_ones[4] = { -1.0, -1.0, -1.0, -1.0 };
418 double loc_side_rhs[4];
428 for (
unsigned int i_loc = 0; i_loc <
el_ds->
lsize(); i_loc++) {
432 unsigned int nsides = ele->n_sides();
436 for (
unsigned int i = 0; i < nsides; i++) {
438 edge_row = edge_rows[i] =
row_4_edge[ele->side(i)->edge_idx()];
439 bcd=ele->side(i)->cond();
442 loc_side_rhs[i] = (ele->centre()[ 2 ] - ele->side(i)->centre()[ 2 ]);
456 loc_side_rhs[i] -= bc_pressure;
480 double val_side = (fe_values.
local_matrix())[i*nsides+i];
481 double val_edge = -1./ (fe_values.
local_matrix())[i*nsides+i];
482 subdomain_diagonal_map.insert( std::make_pair( side_row, val_side ) );
483 subdomain_diagonal_map.insert( std::make_pair( edge_row, val_edge ) );
518 for (
unsigned int i = 0; i < ele->n_neighs_vb; i++) {
521 ngh= ele->neigh_vb[i];
528 switch (ele_higher->dim()) {
543 double value =
data_.
sigma.
value( ele->centre(), ele->element_accessor()) *
552 local_vb[0] = -value; local_vb[1] = value;
553 local_vb[2] = value; local_vb[3] = -value;
559 int ind = tmp_rows[1];
560 double new_val = value;
562 if ( it != subdomain_diagonal_map.end() ) {
563 new_val = new_val + it->second;
565 subdomain_diagonal_map.insert( std::make_pair( ind, new_val ) );
576 tmp_rows[2+i] = tmp_rows[1];
585 ASSERT(ele->n_neighs_vb*ele->n_neighs_vb<1000,
"Too many values in E block.");
587 ele->n_neighs_vb, tmp_rows+2, zeros);
629 vector<int> &dofs,
unsigned int &ele_type,
double &delta, arma::vec &dirichlet) {
633 if (i_ele == ml_it_->size() ) {
639 const Intersection &isect=intersections_[ (*ml_it_)[i_ele] ];
646 for(
unsigned int i_side=0; i_side < ele->
n_sides(); i_side++ ) {
647 dofs[i_side]=darcy_.row_4_edge[ele->
side(i_side)->
edge_idx()];
650 dirichlet.resize(ele->
n_sides());
656 dofs[i_side] = -dofs[i_side];
657 double bc_pressure = darcy_.data_.bc_pressure.value(b_ele.
centre(), b_ele);
658 dirichlet[i_side] = bc_pressure;
669 double delta_i, delta_j;
671 arma::vec dirichlet_i, dirichlet_j;
672 unsigned int ele_type_i, ele_type_j;
677 for(ml_it_ = master_list_.begin(); ml_it_ != master_list_.end(); ++ml_it_) {
679 if (ml_it_->size() == 0)
continue;
695 master_ = intersections_[ml_it_->front()].master_iter();
696 delta_0 = master_->measure();
698 double master_sigma=darcy_.data_.sigma.value( master_->centre(), master_->element_accessor());
701 for(i = 0; i <= ml_it_->size(); ++i) {
704 pressure_diff(i, dofs_i, ele_type_i, delta_i, dirichlet_i);
706 for (j = 0; j <= ml_it_->size(); ++j) {
707 pressure_diff(j, dofs_j, ele_type_j, delta_j, dirichlet_j);
709 double scale = -master_sigma * delta_i * delta_j / delta_0;
710 product = scale * tensor_average[ele_type_i][ele_type_j];
716 arma::vec rhs(dofs_i.size());
718 ls.
set_values( dofs_i, dofs_j, product, rhs, dirichlet_i, dirichlet_j);
729 for(
unsigned int i_side=0; i_side < ele->
n_sides(); i_side++ ) {
730 dofs[shift+i_side] = darcy_.row_4_edge[ele->
side(i_side)->
edge_idx()];
739 dofs[shift + i_side] = -dofs[shift + i_side];
740 double bc_pressure = darcy_.data_.bc_pressure.value(b_ele.
centre(), b_ele);
741 dirichlet[shift + i_side] = bc_pressure;
767 const Element * master = intersec.master_iter();
768 const Element * slave = intersec.slave_iter();
770 add_sides(master, 0, dofs, dirichlet);
771 add_sides(slave, 2, dofs, dirichlet);
798 arma::vec point_Y(1);
800 arma::vec point_2D_Y(intersec.map_to_slave(point_Y));
801 arma::vec point_1D_Y(intersec.map_to_master(point_Y));
803 arma::vec point_X(1);
805 arma::vec point_2D_X(intersec.map_to_slave(point_X));
806 arma::vec point_1D_X(intersec.map_to_master(point_X));
808 arma::mat base_2D(3, 3);
811 base_2D << 1.0 << 0.0 << -2.0 << arma::endr
812 << -1.0 << 2.0 << 2.0 << arma::endr
813 << 1.0 << -2.0 << 0.0 << arma::endr;
815 arma::mat base_1D(2, 2);
818 base_1D << 1.0 << -1.0 << arma::endr
819 << 0.0 << 1.0 << arma::endr;
827 arma::vec difference_in_Y(5);
828 arma::vec difference_in_X(5);
831 difference_in_Y.subvec(0, 2) = -base_2D * point_2D_Y;
832 difference_in_X.subvec(0, 2) = -base_2D * point_2D_X;
834 difference_in_Y.subvec(3, 4) = base_1D * point_1D_Y;
835 difference_in_X.subvec(3, 4) = base_1D * point_1D_X;
842 for (
int i = 0; i < 5; ++i) {
843 for (
int j = 0; j < 5; ++j) {
844 A(i, j) = -master_sigma * intersec.intersection_true_size() *
845 ( difference_in_Y[i] * difference_in_Y[j]
846 + difference_in_Y[i] * difference_in_X[j]/2
847 + difference_in_X[i] * difference_in_Y[j]/2
848 + difference_in_X[i] * difference_in_X[j]
857 ls.
set_values( dofs, dofs, A, rhs, dirichlet, dirichlet);
891 xprintf(
Warn,
"BDDC is using no Schur complements.");
894 xprintf(
Warn,
"Invalid number of Schur Complements. Using 2.");
911 xprintf(
Err,
"Flow123d was not build with BDDCML support.\n");
927 ls->LinSys::set_from_input(in_rec);
935 ASSERT(err == 0,
"Error in ISCreateStride.");
949 ASSERT(err == 0,
"Error in ISCreateStride.");
971 xprintf(
Err,
"Unknown solver type. Internal error.\n");
985 DBGMSG(
"Assembly linear system\n");
987 DBGMSG(
" Data changed\n");
1000 DBGMSG(
" setup time term\n");
1011 xprintf(
PrgErr,
"Planned computation time for steady solver, but data are not changed.\n");
1046 for (
int i_loc = 0; i_loc <
el_ds->
lsize(); i_loc++ ) {
1049 int e_idx = el.
index();
1051 int elDim = el->dim();
1052 elDimMax = std::max( elDimMax, elDim );
1053 elDimMin = std::min( elDimMin, elDim );
1055 isegn.push_back( e_idx );
1063 localDofMap.insert( std::make_pair( side_row, coord ) );
1064 inet.push_back( side_row );
1071 localDofMap.insert( std::make_pair( el_row, coord ) );
1072 inet.push_back( el_row );
1076 Edge *edg=el->side(si)->edge();
1079 int edge_row =
row_4_edge[ el->side(si)->edge_idx() ];
1082 localDofMap.insert( std::make_pair( edge_row, coord ) );
1083 inet.push_back( edge_row );
1088 for (
int i_neigh = 0; i_neigh < el->n_neighs_vb; i_neigh++) {
1089 int edge_row =
row_4_edge[ el->neigh_vb[i_neigh]->edge_idx() ];
1090 arma::vec3 coord = el->neigh_vb[i_neigh]->edge()->side(0)->centre();
1092 localDofMap.insert( std::make_pair( edge_row, coord ) );
1093 inet.push_back( edge_row );
1097 nnet.push_back( nne );
1108 for (
int i = 0; i < 3; i++) {
1109 coef = coef + aniso.at(i,i);
1112 coef = conduct*coef / 3;
1115 "Zero coefficient of hydrodynamic resistance %f . \n ", coef );
1116 element_permeability.push_back( 1. / coef );
1119 int numNodeSub = localDofMap.size();
1125 for ( ; itB != localDofMap.end(); ++itB ) {
1126 isngn[ind] = itB -> first;
1129 for (
int j = 0; j < 3; j++ ) {
1130 xyz[ j*numNodeSub + ind ] = coord[j];
1135 localDofMap.clear();
1142 Global2LocalMap_ global2LocalNodeMap;
1143 for (
unsigned ind = 0; ind < isngn.size(); ++ind ) {
1144 global2LocalNodeMap.insert( std::make_pair( static_cast<unsigned>( isngn[ind] ), ind ) );
1160 for (
int iEle = 0; iEle < isegn.size(); iEle++ ) {
1161 int nne = nnet[ iEle ];
1162 for (
unsigned ien = 0; ien < nne; ien++ ) {
1164 int indGlob = inet[indInet];
1166 Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1167 ASSERT( pos != global2LocalNodeMap.end(),
1168 "Cannot remap node index %d to local indices. \n ", indGlob );
1169 int indLoc =
static_cast<int> ( pos -> second );
1172 inet[ indInet++ ] = indLoc;
1176 int numNodes =
size;
1177 int numDofsInt =
size;
1179 int meshDim = elDimMax;
1182 bddc_ls -> load_mesh( spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1354 int edge_shift[np], el_shift[np], side_shift[np];
1355 unsigned int rows_starts[np];
1363 for (i = 0; i < np; i++) {
1370 rows_starts[i] = shift;
1376 for (i = 0; i < side_n_id; i++) {
1381 for (i = 0; i < el_n_id; i++) {
1387 for (i = 0; i < edge_n_id; i++) {
1393 for (i = np - 1; i > 0; i--)
1394 rows_starts[i] -= rows_starts[i - 1];
1396 rows_ds = boost::make_shared<Distribution>(&(rows_starts[0]), PETSC_COMM_WORLD);
1408 VecCreateSeqWithArray(PETSC_COMM_SELF,1,
size,
solution,
1421 loc_idx[i++] =
row_4_el[ele.index()];
1423 for(
unsigned int i_edg=0; i_edg <
mesh_->
n_edges(); i_edg++) {
1426 ASSERT( i==
size,
"Size of array does not match number of fills.\n");
1428 ISCreateGeneral(PETSC_COMM_SELF,
size, loc_idx, PETSC_COPY_VALUES, &(is_loc));
1432 ISDestroy(&(is_loc));
1477 loc_part =
new int[init_edge_ds.lsize()];
1482 unsigned int i_edg = edg -
mesh_->
edges.begin();
1486 if (init_edge_ds.is_local(i_edg)) {
1488 loc_part[loc_i++] = el_ds->get_proc(
row_4_el[e_idx]);
1491 id_4_old[i_edg] = i_edg;
1502 loc_part =
new int[init_side_ds.lsize()];
1509 if (init_side_ds.is_local(is)) {
1511 loc_part[loc_i++] = el_ds->get_proc(
1550 int side_row, edge_row;
1560 for (
unsigned int i_loc = 0; i_loc < el_ds->lsize(); i_loc++) {
1562 int el_row =
row_4_el[el_4_loc[i_loc]];
1566 unsigned int nsides = el->
n_sides();
1567 for (
unsigned int i = 0; i < nsides; i++) {
1580 for (
unsigned int i_neigh = 0; i_neigh < el->
n_neighs_vb; i_neigh++) {
1588 #endif // HAVE_BDDCML
1601 unsigned int i_edg=0;
1614 const PetscInt *cols;
1620 MatGetOwnershipRange(m, &first, &last);
1621 for (
int row = first; row < last; row++) {
1622 MatGetRow(m, row, &n, &cols, PETSC_NULL);
1623 bool exists_off =
false;
1624 for (
int i = 0; i < n; i++)
1626 n_off++, exists_off =
true;
1631 MatRestoreRow(m, row, &n, &cols, PETSC_NULL);
1692 for (
unsigned int i_loc_el = 0; i_loc_el <
el_ds->
lsize(); i_loc_el++) {
1711 PetscScalar *local_diagonal;
1716 for (
unsigned int i_loc_el = 0; i_loc_el <
el_ds->
lsize(); i_loc_el++) {
1788 for (
unsigned int i_loc_el = 0; i_loc_el <
el_ds->
lsize(); i_loc_el++) {
1794 int edge_row =
row_4_edge[ele->side(i)->edge_idx()];
1820 for (
unsigned int i_loc_el = 0; i_loc_el <
el_ds->
lsize(); i_loc_el++) {
1826 int edge_row =
row_4_edge[ele->side(i)->edge_idx()];
1831 time_->
dt() / ele->n_sides(),ADD_VALUES);
1866 int i_loc, side_row, loc_edge_row, i;
1869 double new_pressure, old_pressure, time_coef;
1871 PetscScalar *loc_prev_sol;
1882 old_pressure = loc_prev_sol[loc_edge_row];
1886 time_coef = - ele->
measure() *
1890 VecSetValue(
schur0->
get_solution(), side_row, time_coef * (new_pressure - old_pressure), ADD_VALUES);
1907 for (
unsigned int i_loc = 0; i_loc <
el_ds->
lsize(); i_loc++) {
1911 values[i] = -1.0 * ele->
measure() *