27 #ifndef OOMPH_BLOCK_PRECONDITION_HEADER
28 #define OOMPH_BLOCK_PRECONDITION_HEADER
33 #include <oomph-lib-config.h>
134 this->
build(0, 0,
false);
149 std::ostringstream err_msg;
150 err_msg <<
"Trying to construct a BlockSelector object with:\n"
151 <<
"replacement_block_pt != 0 and wanted == false"
152 <<
"If you require the block, please set wanted == true.\n";
154 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
173 std::ostringstream warning_msg;
174 warning_msg <<
"Warning: BlockSelector destructor is called but...\n"
175 <<
"replacement_block_pt() is not null.\n"
176 <<
"Please remember to null this via the function\n"
177 <<
"BlockSelector::null_replacement_block_pt()\n"
182 warning_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
196 std::ostringstream err_msg;
197 err_msg <<
"Trying to select block with:\n"
198 <<
"replacement_block_pt != 0 and wanted == false"
199 <<
"If you require the block, please set wanted == true.\n"
203 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
223 std::ostringstream warning_msg;
224 warning_msg <<
"Trying to set Wanted = false, but replacement_block_pt "
226 <<
"Please call null_replacement_block_pt()\n"
227 <<
"(remember to free memory if necessary)\n"
231 warning_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
252 std::ostringstream err_msg;
253 err_msg <<
"Trying to set replacement_block_pt, but Wanted == false.\n"
254 <<
"Please call want_block()\n"
258 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
310 o_stream <<
"Row_index = " << block_selector.
row_index() <<
"\n"
311 <<
"Column_index = " << block_selector.
column_index() <<
"\n"
312 <<
"Wanted = " << block_selector.
wanted() <<
"\n"
313 <<
"Replacement_block_pt = ";
342 std::ostringstream err_msg;
344 <<
"Trying to set replacement_block_pt, but Wanted == false.\n"
345 <<
"Please either not set the replacement_block_pt or call the "
347 <<
"do_not_want_block()\n"
351 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
420 template<
typename MATRIX>
532 std::ostringstream error_msg;
533 error_msg <<
"Matrix is not correct type.";
535 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
675 MATRIX& output_matrix,
676 const bool& ignore_replacement_block =
false)
const
682 if ((i > n_block_types) || (j > n_block_types))
684 std::ostringstream err_msg;
685 err_msg <<
"Requested block(" <<
i <<
"," << j <<
")"
687 <<
"but nblock_types() is " << n_block_types <<
".\n"
688 <<
"Maybe your argument to block_setup(...) is incorrect?"
691 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
780 if (ndofblock_in_row == 1 && ndofblock_in_col == 1)
792 ignore_replacement_block)
797 ignore_replacement_block);
816 ndofblock_in_row, ndofblock_in_col, 0);
829 for (
unsigned row_dofblock = 0; row_dofblock < ndofblock_in_row;
834 const unsigned wanted_dof_row =
837 for (
unsigned col_dofblock = 0; col_dofblock < ndofblock_in_col;
842 const unsigned wanted_dof_col =
846 tmp_blocks_pt(row_dofblock, col_dofblock) =
852 if ((tmp_blocks_pt(row_dofblock, col_dofblock) == 0) ||
853 ignore_replacement_block)
857 new_block[row_dofblock][col_dofblock] = 1;
866 *tmp_blocks_pt(row_dofblock, col_dofblock),
867 ignore_replacement_block);
931 for (
unsigned row_dof_i = 0; row_dof_i < ndofblock_in_row; row_dof_i++)
936 tmp_row_dist_pt[row_dof_i] =
941 tmp_row_dist_pt[row_dof_i] =
950 for (
unsigned col_dof_i = 0; col_dof_i < ndofblock_in_col; col_dof_i++)
955 tmp_col_dist_pt[col_dof_i] =
960 tmp_col_dist_pt[col_dof_i] =
967 tmp_row_dist_pt, tmp_col_dist_pt, tmp_blocks_pt, output_matrix);
970 for (
unsigned row_i = 0; row_i < ndofblock_in_row; row_i++)
972 for (
unsigned col_i = 0; col_i < ndofblock_in_col; col_i++)
974 if (new_block[row_i][col_i])
976 delete tmp_blocks_pt(row_i, col_i);
990 const bool& ignore_replacement_block =
false)
const
992 MATRIX output_matrix;
993 get_block(
i, j, output_matrix, ignore_replacement_block);
994 return output_matrix;
1014 MATRIX* in_matrix_pt,
1015 MATRIX& output_matrix)
1052 MATRIX& output_block,
1053 const bool& ignore_replacement_block =
false)
const;
1139 unsigned const para_selected_block_nrow = selected_block.
nrow();
1140 unsigned const para_selected_block_ncol = selected_block.
ncol();
1144 if (para_selected_block_nrow == 0)
1146 std::ostringstream error_msg;
1147 error_msg <<
"selected_block has nrow = 0.\n";
1149 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1156 if ((para_selected_block_nrow > para_nblocks) ||
1157 (para_selected_block_ncol > para_nblocks))
1159 std::ostringstream error_msg;
1160 error_msg <<
"Trying to concatenate a " << para_selected_block_nrow
1161 <<
" by " << para_selected_block_ncol <<
" block matrix,\n"
1162 <<
"but there are only " << para_nblocks <<
" block types.\n";
1164 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1173 for (
unsigned row_i = 0; row_i < para_selected_block_nrow; row_i++)
1175 const unsigned col_0_row_index = selected_block[row_i][0].row_index();
1177 for (
unsigned col_i = 0; col_i < para_selected_block_ncol; col_i++)
1179 const unsigned para_b_i = selected_block[row_i][col_i].row_index();
1180 const unsigned para_b_j = selected_block[row_i][col_i].column_index();
1182 if (col_0_row_index != para_b_i)
1184 std::ostringstream err_msg;
1185 err_msg <<
"Block index for block(" << row_i <<
"," << col_i <<
") "
1186 <<
"contains block indices (" << para_b_i <<
"," << para_b_j
1188 <<
"But the row index for the first column is "
1189 << col_0_row_index <<
", they must be the same!\n";
1191 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1199 for (
unsigned col_i = 0; col_i < para_selected_block_ncol; col_i++)
1201 const unsigned row_0_col_index =
1202 selected_block[0][col_i].column_index();
1204 for (
unsigned row_i = 0; row_i < para_selected_block_nrow; row_i++)
1206 const unsigned para_b_i = selected_block[row_i][col_i].row_index();
1207 const unsigned para_b_j = selected_block[row_i][col_i].column_index();
1209 if (row_0_col_index != para_b_j)
1211 std::ostringstream err_msg;
1212 err_msg <<
"Block index for block(" << row_i <<
"," << col_i <<
") "
1213 <<
"contains block indices (" << para_b_i <<
"," << para_b_j
1215 <<
"But the col index for the first row is "
1216 << row_0_col_index <<
", they must be the same!\n";
1218 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1232 for (
unsigned row_i = 0; row_i < para_selected_block_nrow; row_i++)
1234 const unsigned para_b_i = selected_block[row_i][0].row_index();
1235 const unsigned para_b_j = selected_block[row_i][0].column_index();
1236 if (para_b_i > para_nblocks)
1238 std::ostringstream err_msg;
1239 err_msg <<
"Block index for block(" << row_i <<
",0) "
1240 <<
"contains block indices (" << para_b_i <<
"," << para_b_j
1242 <<
"But there are only " << para_nblocks
1243 <<
" nblock_types().\n";
1245 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1251 for (
unsigned col_i = 0; col_i < para_selected_block_ncol; col_i++)
1253 const unsigned para_b_i = selected_block[0][col_i].row_index();
1254 const unsigned para_b_j = selected_block[0][col_i].column_index();
1255 if (para_b_j > para_nblocks)
1257 std::ostringstream err_msg;
1258 err_msg <<
"Block index for block(0," << col_i <<
") "
1259 <<
"contains block indices (" << para_b_i <<
"," << para_b_j
1261 <<
"But there are only " << para_nblocks
1262 <<
" nblock_types().\n";
1264 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1270 std::set<unsigned> row_index_set;
1271 for (
unsigned row_i = 0; row_i < para_selected_block_nrow; row_i++)
1273 std::pair<std::set<unsigned>::iterator,
bool> row_index_set_ret;
1275 const unsigned row_i_index = selected_block[row_i][0].row_index();
1277 row_index_set_ret = row_index_set.insert(row_i_index);
1279 if (!row_index_set_ret.second)
1281 std::ostringstream err_msg;
1282 err_msg <<
"The row " << row_i_index <<
" is already inserted.\n";
1284 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1290 std::set<unsigned> col_index_set;
1291 for (
unsigned col_i = 0; col_i < para_selected_block_ncol; col_i++)
1293 std::pair<std::set<unsigned>::iterator,
bool> col_index_set_ret;
1295 const unsigned col_i_index = selected_block[0][col_i].column_index();
1297 col_index_set_ret = col_index_set.insert(col_i_index);
1299 if (!col_index_set_ret.second)
1301 std::ostringstream err_msg;
1302 err_msg <<
"The col " << col_i_index <<
" is already inserted.\n";
1304 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1312 for (
unsigned block_i = 0; block_i < para_selected_block_nrow; block_i++)
1314 for (
unsigned block_j = 0; block_j < para_selected_block_ncol;
1318 selected_block[block_i][block_j].replacement_block_pt();
1320 if (tmp_block_pt != 0)
1322 if (!tmp_block_pt->
built())
1324 std::ostringstream err_msg;
1325 err_msg <<
"The matrix pointed to by block(" << block_i <<
","
1326 << block_j <<
") is not built.\n";
1328 OOMPH_CURRENT_FUNCTION,
1329 OOMPH_EXCEPTION_LOCATION);
1335 const unsigned row_selected_block =
1336 selected_block[block_i][block_j].row_index();
1341 if (*tmp_block_dist_pt != *another_tmp_block_dist_pt)
1343 std::ostringstream err_msg;
1344 err_msg <<
"block_distribution_pt " << row_selected_block <<
"\n"
1345 <<
"does not match the distribution from the block_pt() "
1346 <<
" in selected_block[" << block_i <<
"][" << block_j
1349 OOMPH_CURRENT_FUNCTION,
1350 OOMPH_EXCEPTION_LOCATION);
1365 for (
unsigned block_i = 0; block_i < para_selected_block_nrow; block_i++)
1367 for (
unsigned block_j = 0; block_j < para_selected_block_ncol;
1372 selected_block[block_i][block_j].replacement_block_pt();
1374 if (tmp_block_pt != 0)
1376 const unsigned tmp_block_ncol = tmp_block_pt->
ncol();
1378 const unsigned selected_block_col =
1379 selected_block[block_i][block_j].column_index();
1382 const unsigned another_tmp_block_ncol =
1385 if (tmp_block_ncol != another_tmp_block_ncol)
1387 std::ostringstream err_msg;
1388 err_msg <<
"block_pt in selected_block[" << block_i <<
"]["
1390 <<
"has ncol = " << tmp_block_ncol <<
".\n"
1391 <<
"But the corresponding block_distribution_pt("
1392 << selected_block_col
1393 <<
") has nrow = " << another_tmp_block_ncol <<
").\n";
1395 OOMPH_CURRENT_FUNCTION,
1396 OOMPH_EXCEPTION_LOCATION);
1404 MATRIX output_matrix;
1407 const unsigned nblock_row = selected_block.
nrow();
1408 const unsigned nblock_col = selected_block.
ncol();
1418 for (
unsigned row_i = 0; row_i < nblock_row; row_i++)
1420 const unsigned selected_row_index =
1421 selected_block[row_i][0].row_index();
1424 block_row_index[row_i] = selected_row_index;
1429 for (
unsigned col_i = 0; col_i < nblock_col; col_i++)
1431 const unsigned selected_col_index =
1432 selected_block[0][col_i].column_index();
1435 block_col_index[col_i] = selected_col_index;
1453 output_matrix.build(iter->second);
1461 output_matrix.build(tmp_dist_pt);
1487 for (
unsigned block_i = 0; block_i < nblock_row; block_i++)
1489 for (
unsigned block_j = 0; block_j < nblock_col; block_j++)
1491 const bool block_wanted = selected_block[block_i][block_j].wanted();
1496 selected_block[block_i][block_j].replacement_block_pt();
1498 if (tmp_block_pt == 0)
1500 new_block[block_i][block_j] = 1;
1505 const unsigned tmp_block_i = block_row_index[block_i];
1506 const unsigned tmp_block_j = block_col_index[block_j];
1510 tmp_block_i, tmp_block_j, *block_pt(block_i, block_j));
1514 block_pt(block_i, block_j) = tmp_block_pt;
1522 row_dist_pt, col_dist_pt, block_pt, output_matrix);
1525 for (
unsigned block_i = 0; block_i < nblock_row; block_i++)
1527 for (
unsigned block_j = 0; block_j < nblock_col; block_j++)
1529 if (new_block[block_i][block_j])
1531 delete block_pt(block_i, block_j);
1536 return output_matrix;
1675 std::ostringstream error_msg;
1677 <<
"The Block_to_dof_map_coarse vector is not setup for \n"
1678 <<
"this block preconditioner.\n\n"
1680 <<
"This vector is always set up within block_setup(...).\n"
1681 <<
"If block_setup() is already called, then perhaps there is\n"
1682 <<
"something wrong with your block preconditionable elements.\n"
1685 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1700 std::ostringstream err_msg;
1701 unsigned n =
nmesh();
1704 err_msg <<
"No meshes have been set for this block preconditioner!\n"
1705 <<
"Set one with set_nmesh(...), set_mesh(...)" << std::endl;
1707 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1708 for (
unsigned m = 0; m < n; m++)
1712 err_msg <<
"The mesh pointer to mesh " << m <<
" is null!\n"
1713 <<
"Set a non-null one with set_mesh(...)" << std::endl;
1715 OOMPH_CURRENT_FUNCTION,
1716 OOMPH_EXCEPTION_LOCATION);
1736 std::ostringstream error_msg;
1738 <<
"The Doftype_coarsen_map_coarse vector is not setup for \n"
1739 <<
"this SUBSIDIARY block preconditioner.\n\n"
1741 <<
"For SUBSIDARY block preconditioners at any level, this\n"
1742 <<
"vector is set up in the function \n"
1743 <<
"turn_into_subsidiary_block_preconditioner(...).\n\n"
1745 <<
"Being a SUBSIDIARY block preconditioner means that \n"
1746 <<
"(Master_block_preconditioner_pt == 0) is true.\n"
1747 <<
"The Master_block_preconditioner_pt MUST be set in the "
1749 <<
"turn_into_subsidiary_block_preconditioner(...).\n\n"
1751 <<
"Somewhere, the Master_block_preconditioner_pt pointer is\n"
1752 <<
"set but not by the function\n"
1753 <<
"turn_into_subsidiary_block_preconditioner(...).\n"
1756 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1787 std::ostringstream error_msg;
1788 error_msg <<
"The mesh_pt() function should not be called on a\n"
1789 <<
"subsidiary block preconditioner." << std::endl;
1791 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1800 std::ostringstream error_msg;
1801 error_msg <<
"Mesh pointer " << mesh_i_pt <<
" is null.";
1803 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1834 unsigned block_i = 0;
1855 if (internal_dof_block_number >= 0)
1863 unsigned block_proc =
1868 const unsigned ndof_in_block =
1872 for (
unsigned dof_i = 0; dof_i < ndof_in_block; dof_i++)
1883 internal_dof_block_number)
1904 const unsigned& b)
const
1909 std::ostringstream error_msg;
1910 error_msg <<
"Block distributions are not set up.\n"
1911 <<
"Have you called block_setup(...)?\n"
1914 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1918 std::ostringstream error_msg;
1919 error_msg <<
"You requested the distribution for the block " << b
1922 <<
" block types.\n"
1925 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1938 std::ostringstream error_msg;
1939 error_msg <<
"Block distributions are not set up.\n"
1940 <<
"Have you called block_setup(...)?\n"
1943 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1947 std::ostringstream error_msg;
1948 error_msg <<
"You requested the distribution for the block " << b
1951 <<
" block types.\n"
1954 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1967 std::ostringstream error_msg;
1968 error_msg <<
"You requested the distribution for the dof block " << b
1970 <<
"But there are only " <<
ndof_types() <<
" DOF types.\n"
1973 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1983 std::ostringstream error_msg;
1984 error_msg <<
"Internal block distributions are not set up.\n"
1985 <<
"Have you called block_setup(...)?\n"
1988 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2001 std::ostringstream error_msg;
2002 error_msg <<
"Dof block distributions are not set up.\n"
2003 <<
"Have you called block_setup(...)?\n"
2006 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2042 std::ostringstream err_msg;
2043 err_msg <<
"A subsidiary block preconditioner should not care about\n"
2044 <<
"anything to do with meshes.";
2046 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2096 const unsigned& precision = 8)
const
2100 for (
unsigned i = 0;
i < nblocks;
i++)
2102 for (
unsigned j = 0; j < nblocks; j++)
2110 get_block(
i, j).sparse_indexed_output(filename, precision,
true);
2126 #ifdef OOMPH_HAS_MPI
2146 std::ostringstream error_message;
2147 error_message <<
"This block preconditioner does not have "
2148 <<
"a master preconditioner.";
2150 OOMPH_CURRENT_FUNCTION,
2151 OOMPH_EXCEPTION_LOCATION);
2166 for (
unsigned b = 0; b < nblock; b++)
2178 #ifdef OOMPH_HAS_MPI
2182 for (
unsigned p = 0; p < nc; p++)
2186 for (
unsigned b = 0; b < nr; b++)
2219 for (
unsigned dist_i = 0; dist_i < n_existing_block_dist; dist_i++)
2230 for (
unsigned i = 0;
i < n_existing_block_dist;
i++)
2232 preconditioner_matrix_key[
i] =
i;
2243 if (iter->first != preconditioner_matrix_key)
2245 delete iter->second;
2259 for (
unsigned dof_i = 0; dof_i < ndof_block_dist; dof_i++)
2272 oomph_info <<
"===========================================" << std::endl;
2273 oomph_info <<
"Block Preconditioner Documentation" << std::endl
2284 oomph_info <<
"Master DOF number " << d <<
" : "
2291 oomph_info <<
"Block " << b <<
" DOF types:";
2300 oomph_info <<
"Master block preconditioner distribution:" << std::endl;
2302 oomph_info <<
"Internal preconditioner matrix distribution:" << std::endl;
2305 oomph_info <<
"Preconditioner matrix distribution:" << std::endl;
2309 oomph_info <<
"Internal block " << b <<
" distribution:" << std::endl;
2314 oomph_info <<
"Block " << b <<
" distribution:" << std::endl;
2328 oomph_info <<
"===========================================" << std::endl;
2346 if (
i >= n_dof_types)
2348 std::ostringstream err_msg;
2349 err_msg <<
"Trying to get the most fine grain dof types in dof type "
2350 <<
i <<
",\nbut there are only " << n_dof_types
2351 <<
" number of dof types.\n";
2353 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2366 if (
i >= n_dof_types)
2368 std::ostringstream err_msg;
2369 err_msg <<
"Trying to get the number of most fine grain dof types "
2370 <<
"in dof type " <<
i <<
",\n"
2371 <<
"but there are only " << n_dof_types
2372 <<
" number of dof types.\n";
2374 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2401 const unsigned nblock = block_col_indices.size();
2405 const unsigned col_index = block_col_indices[0];
2416 matvec_prod_pt->
setup(block_pt, iter->second);
2421 for (
unsigned b = 0; b < nblock; b++)
2431 matvec_prod_pt->
setup(block_pt, tmp_dist_pt);
2440 const unsigned& block_col_index)
2527 std::ostringstream err_msg;
2529 <<
"(Internal_nblock_types == 0) is true. \n"
2530 <<
"Did you remember to call the function block_setup(...)?\n\n"
2532 <<
"This variable is always set up within block_setup(...).\n"
2533 <<
"If block_setup() is already called, then perhaps there is\n"
2534 <<
"something wrong with your block preconditionable elements.\n"
2537 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2542 std::ostringstream err_msg;
2544 <<
"The number of internal block types and "
2545 <<
"internal dof types does not match... \n\n"
2546 <<
"Internally, the number of block types and the number of dof "
2547 <<
"types must be the same.\n"
2550 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2573 std::ostringstream error_msg;
2575 <<
"(Internal_ndof_types == 0) is true.\n"
2576 <<
"This means that the Master_block_preconditioner_pt pointer is\n"
2577 <<
"set but possibly not by the function\n"
2578 <<
"turn_into_subsidiary_block_preconditioner(...).\n\n"
2580 <<
"This goes against the block preconditioning framework "
2582 <<
"Many machinery relies on the look up lists set up by the \n"
2583 <<
"function turn_into_subsidiary_block_preconditioner(...) \n"
2584 <<
"between the parent and child block preconditioners.\n"
2587 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2597 for (
unsigned i = 0;
i <
nmesh();
i++)
2679 MATRIX& output_block)
const;
2712 if (internal_dof_block_number >= 0)
2720 internal_dof_block_number)
2735 const unsigned& b)
const
2740 std::ostringstream error_msg;
2741 error_msg <<
"Internal block distributions are not set up.\n"
2742 <<
"Have you called block_setup(...)?\n"
2745 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2749 std::ostringstream error_msg;
2750 error_msg <<
"You requested the distribution for the internal block "
2753 <<
" block types.\n"
2756 error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2774 const unsigned max_block_number =
2775 *std::max_element(block_vec_number.begin(), block_vec_number.end());
2778 if (max_block_number >= nblocks)
2780 std::ostringstream err_msg;
2781 err_msg <<
"Cannot insert into Auxiliary_block_distribution_pt\n"
2782 <<
"because " << max_block_number <<
" is equal to or \n"
2783 <<
"greater than " << nblocks <<
".\n";
2785 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2800 std::ostringstream err_msg;
2801 err_msg <<
"Cannot insert into Auxiliary_block_distribution_pt\n"
2802 <<
"because the first in the pair already exists.\n";
2804 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2809 std::make_pair(block_vec_number, dist_pt));
2816 const MATRIX* block_matrix_pt)
const;
2826 template<
typename myType>
2829 const bool sorted =
false)
const
2834 std::lower_bound(vec.begin(), vec.end(), val);
2836 return (low == vec.end() || *low != val) ? -1 : (low - vec.begin());
2840 int pos = std::find(vec.begin(), vec.end(), val) - vec.begin();
2841 return (pos <
int(vec.size()) && pos >= 0) ? pos : -1;
2868 const bool& allow_multiple_element_type_in_mesh =
false)
2874 std::ostringstream err_msg;
2875 err_msg <<
"The mesh pointer has space for " <<
nmesh() <<
" meshes.\n"
2876 <<
"Cannot store a mesh at entry " <<
i <<
"\n"
2877 <<
"Has set_nmesh(...) been called?";
2879 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2885 std::ostringstream err_msg;
2886 err_msg <<
"Tried to set the " <<
i
2887 <<
"-th mesh pointer, but it is null.";
2889 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2898 unsigned(allow_multiple_element_type_in_mesh);
2912 const unsigned& block_j,
2919 std::ostringstream err_msg;
2920 err_msg <<
"nblock_types() is 0, has block_setup(...) been called?\n";
2922 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2927 unsigned para_ndof_types = this->
ndof_types();
2929 if ((block_i >= para_ndof_types) || (block_j >= para_ndof_types))
2931 std::ostringstream err_msg;
2932 err_msg <<
"Replacement dof block (" << block_i <<
"," << block_j
2933 <<
") is outside of range:\n"
2936 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2957 std::ostringstream err_msg;
2958 err_msg <<
"Replacing block(" << block_i <<
"," << block_i <<
")\n"
2959 <<
" but the pointer is NULL." << std::endl;
2961 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2967 std::ostringstream err_msg;
2968 err_msg <<
"Replacement block(" << block_i <<
"," << block_i <<
")"
2969 <<
" is not built." << std::endl;
2971 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2979 const unsigned para_dof_block_i = block_i;
2984 std::ostringstream err_msg;
2985 err_msg <<
"The distribution of the replacement dof_block_pt\n"
2986 <<
"is different from the Dof_block_distribution_pt["
2987 << para_dof_block_i <<
"].\n";
2989 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2994 const unsigned para_dof_block_j = block_j;
2996 unsigned para_required_ncol =
2998 if (para_replacement_block_ncol != para_required_ncol)
3000 std::ostringstream err_msg;
3001 err_msg <<
"Replacement dof block has ncol = "
3002 << para_replacement_block_ncol <<
".\n"
3003 <<
"But required ncol is " << para_required_ncol <<
".\n";
3005 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
3027 #ifdef OOMPH_HAS_MPI
3029 for (
unsigned i = 0, n =
nmesh();
i < n;
i++)
3031 if (
mesh_pt(
i)->is_mesh_distributed())
3049 #ifdef OOMPH_HAS_MPI
3053 if (i_dof >=
first_row && i_dof <= last_row)
3069 unsigned my_rank =
comm_pt()->my_rank();
3070 std::ostringstream error_message;
3071 error_message <<
"Proc " << my_rank
3072 <<
": Requested internal_dof_number(...) for global DOF "
3074 <<
"cannot be found.\n";
3076 OOMPH_CURRENT_FUNCTION,
3077 OOMPH_EXCEPTION_LOCATION);
3097 return static_cast<int>(
i);
3106 "Never get here\n", OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
3117 #ifdef OOMPH_HAS_MPI
3121 if (i_dof >=
first_row && i_dof <= last_row)
3137 std::ostringstream error_message;
3138 error_message <<
"Requested internal_index_in_dof(...) for global DOF "
3140 <<
"cannot be found.\n";
3142 OOMPH_CURRENT_FUNCTION,
3143 OOMPH_EXCEPTION_LOCATION);
3156 "Never get here\n", OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
3169 if (b >= i_nblock_types)
3171 std::ostringstream err_msg;
3172 err_msg <<
"Trying to get internal block dimension for \n"
3173 <<
"internal block " << b <<
".\n"
3174 <<
"But there are only " << i_nblock_types
3175 <<
" internal dof types.\n";
3177 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
3191 if (
i >= i_n_dof_types)
3193 std::ostringstream err_msg;
3194 err_msg <<
"Trying to get internal dof block dimension for \n"
3195 <<
"internal dof block " <<
i <<
".\n"
3196 <<
"But there are only " << i_n_dof_types
3197 <<
" internal dof types.\n";
3199 err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
3422 #ifdef OOMPH_HAS_MPI
3477 #ifdef OOMPH_HAS_MPI
Block Preconditioner base class. The block structure of the overall problem is determined from the Me...
const LinearAlgebraDistribution * master_distribution_pt() const
Access function to the distribution of the master preconditioner. If this preconditioner does not hav...
Vector< Vector< unsigned > > Doftype_coarsen_map_coarse
Mapping for dof types within THIS precondition. This is usually passed down from the parent precondit...
Vector< unsigned > Nrows_to_recv_for_get_ordered
The number of preconditioner rows to be received from processor p for get_block_ordered_....
void internal_return_block_vector(const unsigned &n, const DoubleVector &b, DoubleVector &v) const
Takes the n-th block ordered vector, b, and copies its entries to the appropriate entries in the natu...
BlockPreconditioner(const BlockPreconditioner &)=delete
Broken copy constructor.
Vector< unsigned > Index_in_dof_block_sparse
Vector to store the mapping from the global DOF number to the index (row/column number) within its bl...
unsigned Internal_nblock_types
Number of different block types in this preconditioner. Note that this information is maintained if u...
Vector< unsigned > Index_in_dof_block_dense
This was uncommented Presumably a non-distributed analogue of Index_in_dof_block_sparse.
bool block_output_on() const
Test if output of blocks is on or not.
void turn_on_recursive_debug_flag()
Toggles on the recursive debug flag. The change goes up the block preconditioning hierarchy.
unsigned nmesh() const
Return the number of meshes in Mesh_pt.
const LinearAlgebraDistribution * internal_block_distribution_pt(const unsigned &b) const
Access function to the internal block distributions.
void return_block_vector(const unsigned &n, const DoubleVector &b, DoubleVector &v) const
Takes the n-th block ordered vector, b, and copies its entries to the appropriate entries in the natu...
unsigned internal_index_in_dof(const unsigned &i_dof) const
Return the row/column number of global unknown i_dof within it's block.
LinearAlgebraDistribution * dof_block_distribution_pt(const unsigned &b)
Access function to the dof-level block distributions.
DenseMatrix< unsigned > Nrows_to_send_for_get_block
The number of global rows to be sent of block b to processor p (matrix indexed [b][p]).
Vector< LinearAlgebraDistribution * > Dof_block_distribution_pt
Storage for the default distribution for each dof block at this level.
void turn_on_debug_flag()
Toggles on the debug flag.
Vector< Vector< unsigned > > Block_to_dof_map_fine
Mapping for the block types to the most fine grain dof types.
void document()
debugging method to document the setup. Should only be called after block_setup(.....
void internal_return_block_vectors(const Vector< unsigned > &block_vec_number, const Vector< DoubleVector > &s, DoubleVector &v) const
A helper function, takes the vector of block vectors, s, and copies its entries into the naturally or...
Vector< unsigned > Doftype_in_master_preconditioner_coarse
The map between the dof types in this preconditioner and the master preconditioner....
Vector< unsigned > Allow_multiple_element_type_in_mesh
Vector of unsigned to indicate which meshes contain multiple element types.
void get_dof_level_block(const unsigned &i, const unsigned &j, MATRIX &output_block, const bool &ignore_replacement_block=false) const
Gets dof-level block (i,j). If Replacement_dof_block_pt(i,j) is not null, then the replacement block ...
Vector< unsigned > Dof_number_dense
Vector to store the mapping from the global DOF number to its block. Empty if this preconditioner has...
DenseMatrix< int * > Rows_to_recv_for_get_block
The block rows to be received from processor p for block b (matrix indexed [b][p]).
Vector< LinearAlgebraDistribution * > Internal_block_distribution_pt
Storage for the default distribution for each internal block.
unsigned Internal_ndof_types
Number of different DOF types in this preconditioner. Note that this information is maintained if use...
const LinearAlgebraDistribution * internal_preconditioner_matrix_distribution_pt() const
access function to the internal preconditioner matrix distribution pt. preconditioner_matrix_distribu...
Vector< Vector< unsigned > > doftype_coarsen_map_fine() const
Access function for the Doftype_coarsen_map_fine variable.
const Mesh * mesh_pt(const unsigned &i) const
Access to i-th mesh (of the various meshes that contain block preconditionable elements of the same n...
unsigned master_nrow() const
Return the number of dofs (number of rows or columns) in the overall problem. The prefix "master_" is...
Vector< unsigned > get_fine_grain_dof_types_in(const unsigned &i) const
Returns the most fine grain dof types in a (possibly coarsened) dof type.
Vector< unsigned > Ndof_types_in_mesh
Storage for number of types of degree of freedom of the elements in each mesh.
unsigned ndof_types_in_mesh(const unsigned &i) const
Return the number of DOF types in mesh i. WARNING: This should only be used by the upper-most master ...
void disable_block_output_to_files()
Turn off output of blocks (by clearing the basefilename string).
void turn_off_debug_flag()
Toggles off the debug flag.
Vector< int * > Rows_to_send_for_get_ordered
The global rows to be sent to processor p for get_block_ordered_... type methods.
unsigned Nrow
Number of DOFs (# of rows or columns in the matrix) in this preconditioner. Note that this informatio...
MATRIX get_block(const unsigned &i, const unsigned &j, const bool &ignore_replacement_block=false) const
Return block (i,j). If the optional argument ignore_replacement_block is true, then any blocks in Rep...
void block_matrix_test(const unsigned &i, const unsigned &j, const MATRIX *block_matrix_pt) const
Private helper function to check that every element in the block matrix (i,j) matches the correspondi...
unsigned internal_ndof_types() const
Return the number of internal dof types. This is the number of most fine grain dof types....
unsigned nfine_grain_dof_types_in(const unsigned &i) const
Access function for the number of most fine grain dof types in a (possibly coarsened) dof type.
MATRIX get_concatenated_block(const VectorMatrix< BlockSelector > &selected_block)
Returns a concatenation of the block matrices specified by the argument selected_block....
Vector< unsigned > Dof_dimension
Vector containing the size of each block, i.e. the number of global DOFs associated with it....
MapMatrix< unsigned, CRDoubleMatrix * > replacement_dof_block_pt() const
Access function to the replaced dof-level blocks.
LinearAlgebraDistribution * Internal_preconditioner_matrix_distribution_pt
The distribution of the (internal) preconditioner matrix. This is formed by concatenating the distrib...
Vector< LinearAlgebraDistribution * > Block_distribution_pt
The distribution for the blocks.
Vector< unsigned > Dof_number_sparse
Vector to store the mapping from the global DOF number to its block (empty if this preconditioner has...
Vector< Vector< unsigned > > Block_number_to_dof_number_lookup
Vector of vectors to store the mapping from block number to the DOF number (each element could be a v...
void setup_matrix_vector_product(MatrixVectorProduct *matvec_prod_pt, CRDoubleMatrix *block_pt, const unsigned &block_col_index)
Setup matrix vector product. This is simply a wrapper around the other setup_matrix_vector_product fu...
void internal_return_block_ordered_preconditioner_vector(const DoubleVector &w, DoubleVector &v) const
Takes the block ordered vector, w, and reorders it in the natural order. Reordered vector is returned...
Vector< int * > Rows_to_recv_for_get_ordered
The preconditioner rows to be received from processor p for get_block_ordered_... type methods.
void get_blocks(DenseMatrix< bool > &required_blocks, DenseMatrix< MATRIX * > &block_matrix_pt) const
Get all the block matrices required by the block preconditioner. Takes a pointer to a matrix of bools...
void return_concatenated_block_vector(const Vector< unsigned > &block_vec_number, const DoubleVector &b, DoubleVector &v) const
Takes concatenated block ordered vector, b, and copies its entries to the appropriate entries in the ...
static bool Run_block_matrix_test
Static boolean to allow block_matrix_test(...) to be run. Defaults to false.
bool any_mesh_distributed() const
Check if any of the meshes are distributed. This is equivalent to problem.distributed() and is used a...
void return_block_vectors(const Vector< unsigned > &block_vec_number, const Vector< DoubleVector > &s, DoubleVector &v) const
Takes the vector of block vectors, s, and copies its entries into the naturally ordered vector,...
unsigned internal_master_dof_number(const unsigned &b) const
Takes the block number within this preconditioner and returns the corresponding block number in the m...
void get_block(const unsigned &i, const unsigned &j, MATRIX &output_matrix, const bool &ignore_replacement_block=false) const
Put block (i,j) into output_matrix. This block accounts for any coarsening of dof types and any repla...
void set_master_matrix_pt(MATRIX *in_matrix_pt)
Set the matrix_pt in the upper-most master preconditioner.
BlockPreconditioner()
Constructor.
const LinearAlgebraDistribution * preconditioner_matrix_distribution_pt() const
Access function to the preconditioner matrix distribution pointer. This is the concatenation of the b...
int index_in_block(const unsigned &i_dof) const
Given a global dof number, returns the index in the block it belongs to. This is the overall index,...
BlockPreconditioner< MATRIX > * Master_block_preconditioner_pt
If the block preconditioner is acting a subsidiary block preconditioner then a pointer to the master ...
const LinearAlgebraDistribution * block_distribution_pt(const unsigned &b) const
Access function to the block distributions (const version).
void clear_block_preconditioner_base()
Clears all BlockPreconditioner data. Called by the destructor and the block_setup(....
virtual ~BlockPreconditioner()
Destructor.
Vector< Vector< unsigned > > Global_index
Vectors of vectors for the mapping from block number and block row to global row number....
unsigned internal_nblock_types() const
Return the number internal blocks. This should be the same as the number of internal dof types....
void get_block_vector(const unsigned &n, const DoubleVector &v, DoubleVector &b) const
Takes the naturally ordered vector, v and returns the n-th block vector, b. Here n is the block numbe...
bool is_master_block_preconditioner() const
Return true if this preconditioner is the master block preconditioner.
unsigned nblock_types() const
Return the number of block types.
Vector< unsigned > Ndof_in_block
Number of types of degree of freedom associated with each block.
Vector< const Mesh * > Mesh_pt
Vector of pointers to the meshes containing the elements used in the block preconditioner....
void set_replacement_dof_block(const unsigned &block_i, const unsigned &block_j, CRDoubleMatrix *replacement_dof_block_pt)
Set replacement dof-level blocks. Only dof-level blocks can be set. This is important due to how the ...
int block_number(const unsigned &i_dof) const
Return the block number corresponding to a global index i_dof.
void turn_off_recursive_debug_flag()
Toggles off the recursive debug flag. The change goes up the block preconditioning hierarchy.
void operator=(const BlockPreconditioner &)=delete
Broken assignment operator.
std::string Output_base_filename
String giving the base of the files to write block data into. If empty then do not output blocks....
void return_block_ordered_preconditioner_vector(const DoubleVector &w, DoubleVector &v) const
Takes the block ordered vector, w, and reorders it in natural order. Reordered vector is returned in ...
MATRIX * matrix_pt() const
Access function to matrix_pt. If this is the master then cast the matrix pointer to MATRIX*,...
DenseMatrix< unsigned > Nrows_to_recv_for_get_block
The number of block rows to be received from processor p for block b (matrix indexed [b][p]).
bool is_subsidiary_block_preconditioner() const
Return true if this preconditioner is a subsidiary preconditioner.
void insert_auxiliary_block_distribution(const Vector< unsigned > &block_vec_number, LinearAlgebraDistribution *dist_pt)
insert a Vector<unsigned> and LinearAlgebraDistribution* pair into Auxiliary_block_distribution_pt....
void get_block_other_matrix(const unsigned &i, const unsigned &j, MATRIX *in_matrix_pt, MATRIX &output_matrix)
Get a block from a different matrix using the blocking scheme that has already been set up.
Vector< Vector< unsigned > > Doftype_coarsen_map_fine
Mapping the dof types within this preconditioner. The values in here refers to the most grain dof typ...
void post_block_matrix_assembly_partial_clear()
A helper method to reduce the memory requirements of block preconditioners. Once the methods get_bloc...
std::map< Vector< unsigned >, LinearAlgebraDistribution * > Auxiliary_block_distribution_pt
Stores any block-level distributions / concatenation of block-level distributions required....
void set_block_output_to_files(const std::string &basefilename)
Set the base part of the filename to output blocks to. If it is set then all blocks will be output at...
bool Debug_flag
Debugging variable. Set true or false via the access functions turn_on_debug_flag(....
unsigned ndof_types() const
Return the total number of DOF types.
BlockPreconditioner< MATRIX > * master_block_preconditioner_pt() const
Access function to the master block preconditioner pt.
void turn_into_subsidiary_block_preconditioner(BlockPreconditioner< MATRIX > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse)
Function to turn this preconditioner into a subsidiary preconditioner that operates within a bigger "...
int internal_dof_number(const unsigned &i_dof) const
Return the number of the block associated with global unknown i_dof. If this preconditioner is a subs...
void internal_get_block(const unsigned &i, const unsigned &j, MATRIX &output_block) const
Gets block (i,j) from the matrix pointed to by Matrix_pt and returns it in output_block....
unsigned internal_dof_block_dimension(const unsigned &i) const
Return the size of the dof "block" i, i.e. how many degrees of freedom are associated with it....
DenseMatrix< int * > Rows_to_send_for_get_block
The global rows to be sent of block b to processor p (matrix indexed [b][p]).
unsigned internal_block_dimension(const unsigned &b) const
Return the number of degrees of freedom in block b. Note that if this preconditioner acts as a subsid...
int get_index_of_value(const Vector< myType > &vec, const myType val, const bool sorted=false) const
Get the index of first occurrence of value in a vector. If the element does not exist,...
LinearAlgebraDistribution * block_distribution_pt(const unsigned b)
Access function to the block distributions (non-const version).
void internal_get_block_vectors(const Vector< unsigned > &block_vec_number, const DoubleVector &v, Vector< DoubleVector > &s) const
Takes the naturally ordered vector and rearranges it into a vector of sub vectors corresponding to th...
MapMatrix< unsigned, CRDoubleMatrix * > Replacement_dof_block_pt
The replacement dof-level blocks.
void get_block_vectors(const Vector< unsigned > &block_vec_number, const DoubleVector &v, Vector< DoubleVector > &s) const
Takes the naturally ordered vector and rearranges it into a vector of sub vectors corresponding to th...
LinearAlgebraDistribution * Preconditioner_matrix_distribution_pt
The distribution of the preconditioner matrix. This is the concatenation of the block distribution.
void setup_matrix_vector_product(MatrixVectorProduct *matvec_prod_pt, CRDoubleMatrix *block_pt, const Vector< unsigned > &block_col_indices)
Setup a matrix vector product. matvec_prod_pt is a pointer to the MatrixVectorProduct,...
void output_blocks_to_files(const std::string &basefilename, const unsigned &precision=8) const
Output all blocks to numbered files. Called at the end of get blocks if an output filename has been s...
Vector< Vector< unsigned > > Block_to_dof_map_coarse
Mapping for block types to dof types. These are the dof types the writer of the preconditioner expect...
void set_nmesh(const unsigned &n)
Specify the number of meshes required by this block preconditioner. Note: elements in different meshe...
Vector< unsigned > Nrows_to_send_for_get_ordered
The number global rows to be sent to processor p for get_block_ordered_... type methods.
Vector< unsigned > Dof_number_to_block_number_lookup
Vector to the mapping from DOF number to block number.
void internal_get_block_vector(const unsigned &n, const DoubleVector &v, DoubleVector &b) const
A helper function, takes the naturally ordered vector, v, and extracts the n-th block vector,...
bool Recursive_debug_flag
Debugging variable. Set true or false via the access functions turn_on_recursive_debug_flag(....
Vector< unsigned > Global_index_sparse
For global indices outside of the range this->first_row() to this->first_row()+this->nrow_local(),...
Vector< unsigned > Doftype_in_master_preconditioner_fine
The map between the dof types in this preconditioner and the master preconditioner....
int internal_block_number(const unsigned &i_dof) const
Return the block number corresponding to a global index i_dof. This returns the block number correspo...
virtual void block_setup()
Determine the size of the matrix blocks and setup the lookup schemes relating the global degrees of f...
void get_concatenated_block_vector(const Vector< unsigned > &block_vec_number, const DoubleVector &v, DoubleVector &b)
Takes the naturally ordered vector and extracts the blocks indicated by the block number (the values)...
void set_mesh(const unsigned &i, const Mesh *const mesh_pt, const bool &allow_multiple_element_type_in_mesh=false)
Set the i-th mesh for this block preconditioner. Note: The method set_nmesh(...) must be called befor...
void internal_get_block_ordered_preconditioner_vector(const DoubleVector &v, DoubleVector &w) const
Given the naturally ordered vector, v, return the vector rearranged in block order in w....
void get_block_ordered_preconditioner_vector(const DoubleVector &v, DoubleVector &w)
Given the naturally ordered vector, v, return the vector rearranged in block order in w....
int internal_index_in_block(const unsigned &i_dof) const
Return the index in the block corresponding to a global block number i_dof. The index returned corres...
Data structure to store information about a certain "block" or sub-matrix from the overall matrix in ...
void null_replacement_block_pt()
Set Replacement_block_pt to null.
void do_not_want_block()
Indicate that we do not want the block (set Wanted to false).
void want_block()
Indicate that we require the block (set Wanted to true).
const unsigned & column_index() const
returns the column index.
virtual ~BlockSelector()
Default destructor.
void build(const unsigned &row_index, const unsigned &column_index, const bool &wanted, CRDoubleMatrix *replacement_block_pt=0)
Build function, sets the Row_index, Column_index and Wanted variables. the Replacement_block_pt is on...
void select_block(const unsigned &row_index, const unsigned &column_index, const bool &wanted, CRDoubleMatrix *replacement_block_pt=0)
Select a block.
BlockSelector(const unsigned &row_index, const unsigned &column_index, const bool &wanted, CRDoubleMatrix *replacement_block_pt=0)
Constructor, takes the row and column indices and a boolean indicating if the block is required or no...
friend std::ostream & operator<<(std::ostream &o_stream, const BlockSelector &block_selector)
Output function, outputs the Row_index, Column_index, Wanted and the address of the Replacement_block...
BlockSelector()
Default constructor, initialise block index i, j to 0 and bool to false.
const unsigned & row_index() const
returns the row index.
void set_column_index(const unsigned &column_index)
Set the column index.
const bool & wanted() const
returns whether the block is wanted or not.
void set_row_index(const unsigned &row_index)
Set the row index.
CRDoubleMatrix * Replacement_block_pt
Pointer to the block.
CRDoubleMatrix * replacement_block_pt() const
Returns Replacement_block_pt.
bool Wanted
Bool to indicate if we require this block.
unsigned Column_index
Column index of the block.
unsigned Row_index
Row index of the block.
void set_replacement_block_pt(CRDoubleMatrix *replacement_block_pt)
set Replacement_block_pt.
A class for compressed row matrices. This is a distributable object.
unsigned long ncol() const
Return the number of columns of the matrix.
bool built() const
access function to the Built flag - indicates whether the matrix has been build - i....
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
unsigned long nrow() const
Return the number of rows of the matrix.
unsigned long ncol() const
Return the number of columns of the matrix.
void resize(const unsigned long &n)
Resize to a square nxn matrix; any values already present will be transfered.
void clear_distribution()
clear the distribution of this distributable linear algebra object
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
unsigned nrow_local() const
access function for the num of local rows on this processor.
unsigned first_row() const
access function for the first row on this processor
A vector in the mathematical sense, initially developed for linear algebra type applications....
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
unsigned first_row() const
access function for the first row on this processor. If not distributed then this is just zero.
unsigned nrow() const
access function to the number of global rows.
unsigned nrow_local() const
access function for the num of local rows on this processor. If no MPI then Nrow is returned.
unsigned rank_of_global_row(const unsigned i) const
return the processor rank of the global row number i
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Matrix vector product helper class - primarily a wrapper to Trilinos's Epetra matrix vector product m...
void setup(CRDoubleMatrix *matrix_pt, const LinearAlgebraDistribution *col_dist_pt=0)
Setup the matrix vector product operator. WARNING: This class is wrapper to Trilinos Epetra matrix ve...
unsigned ndof_types() const
Return number of dof types in mesh.
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
virtual DoubleMatrixBase * matrix_pt() const
Get function for matrix pointer.
virtual const OomphCommunicator * comm_pt() const
Get function for comm pointer.
virtual void set_matrix_pt(DoubleMatrixBase *matrix_pt)
Set the matrix pointer.
VectorMatrix is a generalised, STL-map-based, matrix based on a Vector of Vectors.
const unsigned ncol() const
return the number of columns. This is the size of the first inner vectors, or returns 0 if the outer ...
const unsigned nrow() const
returns the number of rows. This is the outer Vector size.
void deep_copy(const CRDoubleMatrix *const in_matrix_pt, CRDoubleMatrix &out_matrix)
Create a deep copy of the matrix pointed to by in_matrix_pt.
void concatenate_without_communication(const Vector< LinearAlgebraDistribution * > &row_distribution_pt, const Vector< LinearAlgebraDistribution * > &col_distribution_pt, const DenseMatrix< CRDoubleMatrix * > &matrix_pt, CRDoubleMatrix &result_matrix)
Concatenate CRDoubleMatrix matrices.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void concatenate(const Vector< LinearAlgebraDistribution * > &in_distribution_pt, LinearAlgebraDistribution &out_distribution)
Takes a vector of LinearAlgebraDistribution objects and concatenates them such that the nrow_local of...
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types).
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...