17 #ifndef dealii__matrix_free_h 18 #define dealii__matrix_free_h 20 #include <deal.II/base/exceptions.h> 21 #include <deal.II/base/parallel.h> 22 #include <deal.II/base/quadrature.h> 23 #include <deal.II/base/vectorization.h> 24 #include <deal.II/base/template_constraints.h> 25 #include <deal.II/fe/fe.h> 26 #include <deal.II/fe/mapping.h> 27 #include <deal.II/fe/mapping_q1.h> 28 #include <deal.II/lac/vector.h> 29 #include <deal.II/lac/parallel_vector.h> 30 #include <deal.II/lac/block_vector_base.h> 31 #include <deal.II/lac/constraint_matrix.h> 32 #include <deal.II/dofs/dof_handler.h> 33 #include <deal.II/hp/dof_handler.h> 34 #include <deal.II/hp/q_collection.h> 35 #include <deal.II/matrix_free/helper_functions.h> 36 #include <deal.II/matrix_free/shape_info.h> 37 #include <deal.II/matrix_free/dof_info.h> 38 #include <deal.II/matrix_free/mapping_info.h> 40 #ifdef DEAL_II_WITH_THREADS 42 #include <tbb/task_scheduler_init.h> 43 #include <tbb/parallel_for.h> 44 #include <tbb/blocked_range.h> 52 DEAL_II_NAMESPACE_OPEN
103 template <
int dim,
typename Number=
double>
295 template <
typename DoFHandlerType,
typename QuadratureType>
297 const DoFHandlerType &dof_handler,
300 const QuadratureType &quad,
307 template <
typename DoFHandlerType,
typename QuadratureType>
309 const DoFHandlerType &dof_handler,
311 const QuadratureType &quad,
318 template <
typename DoFHandlerType,
typename QuadratureType>
319 void reinit (
const DoFHandlerType &dof_handler,
321 const QuadratureType &quad,
350 template <
typename DoFHandlerType,
typename QuadratureType>
352 const std::vector<const DoFHandlerType *> &dof_handler,
353 const std::vector<const ConstraintMatrix *> &constraint,
354 const std::vector<IndexSet> &locally_owned_set,
355 const std::vector<QuadratureType> &quad,
363 template <
typename DoFHandlerType,
typename QuadratureType>
365 const std::vector<const DoFHandlerType *> &dof_handler,
366 const std::vector<const ConstraintMatrix *> &constraint,
367 const std::vector<QuadratureType> &quad,
374 template <
typename DoFHandlerType,
typename QuadratureType>
375 void reinit (
const std::vector<const DoFHandlerType *> &dof_handler,
376 const std::vector<const ConstraintMatrix *> &constraint,
377 const std::vector<QuadratureType> &quad,
387 template <
typename DoFHandlerType,
typename QuadratureType>
389 const std::vector<const DoFHandlerType *> &dof_handler,
390 const std::vector<const ConstraintMatrix *> &constraint,
391 const QuadratureType &quad,
398 template <
typename DoFHandlerType,
typename QuadratureType>
399 void reinit (
const std::vector<const DoFHandlerType *> &dof_handler,
400 const std::vector<const ConstraintMatrix *> &constraint,
401 const QuadratureType &quad,
436 template <
typename OutVector,
typename InVector>
440 const std::pair<
unsigned int,
441 unsigned int> &)> &cell_operation,
443 const InVector &src)
const;
454 template <
typename CLASS,
typename OutVector,
typename InVector>
458 const std::pair<
unsigned int,
459 unsigned int> &)
const,
460 const CLASS *owning_class,
462 const InVector &src)
const;
467 template <
typename CLASS,
typename OutVector,
typename InVector>
471 const std::pair<
unsigned int,
475 const InVector &src)
const;
484 std::pair<unsigned int,unsigned int>
486 const unsigned int fe_degree,
487 const unsigned int vector_component = 0)
const;
495 std::pair<unsigned int,unsigned int>
497 const unsigned int fe_index,
498 const unsigned int vector_component = 0)
const;
514 template <
typename VectorType>
516 const unsigned int vector_component=0)
const;
526 template <
typename Number2>
528 const unsigned int vector_component=0)
const;
540 const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &
564 const std::vector<unsigned int> &
571 void renumber_dofs (std::vector<types::global_dof_index> &renumbering,
572 const unsigned int vector_component = 0);
626 const unsigned int vector_number,
627 const unsigned int fe_component = 0)
const;
642 const unsigned int vector_number,
643 const unsigned int fe_component = 0)
const;
675 const unsigned int hp_active_fe_index = 0)
const;
682 const unsigned int hp_active_fe_index = 0)
const;
690 const unsigned int hp_active_fe_index = 0)
const;
698 const unsigned int hp_active_fe_index = 0)
const;
705 const unsigned int hp_active_fe_index = 0)
const;
712 const unsigned int hp_active_fe_index = 0)
const;
736 template <
typename StreamType>
743 void print (std::ostream &out)
const;
767 get_mapping_info ()
const;
772 const internal::MatrixFreeFunctions::DoFInfo &
773 get_dof_info (
const unsigned int fe_component = 0)
const;
800 const unsigned int quad_index = 0,
801 const unsigned int hp_active_fe_index = 0,
802 const unsigned int hp_active_quad_index = 0)
const;
814 const std::vector<const ConstraintMatrix *> &constraint,
815 const std::vector<IndexSet> &locally_owned_set,
824 const std::vector<const ConstraintMatrix *> &constraint,
825 const std::vector<IndexSet> &locally_owned_set,
837 const std::vector<IndexSet> &locally_owned_set);
843 const unsigned int level);
849 const unsigned int level);
859 std::vector<SmartPointer<const DoFHandler<dim> > > dof_handler;
860 std::vector<SmartPointer<const hp::DoFHandler<dim> > > hp_dof_handler;
861 enum ActiveDoFHandler { usual,
hp } active_dof_handler;
862 unsigned int n_dof_handlers;
875 std::vector<internal::MatrixFreeFunctions::DoFInfo>
dof_info;
939 template <
int dim,
typename Number>
940 template <
typename VectorType>
944 const unsigned int comp)
const 947 vec.reinit(
dof_info[comp].vector_partitioner->size());
952 template <
int dim,
typename Number>
953 template <
typename Number2>
957 const unsigned int comp)
const 965 template <
int dim,
typename Number>
967 const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> &
971 return dof_info[comp].vector_partitioner;
976 template <
int dim,
typename Number>
978 const std::vector<unsigned int> &
982 return dof_info[comp].constrained_dofs;
987 template <
int dim,
typename Number>
998 template <
int dim,
typename Number>
1008 template <
int dim,
typename Number>
1018 template <
int dim,
typename Number>
1028 template <
int dim,
typename Number>
1038 template <
int dim,
typename Number>
1048 template <
int dim,
typename Number>
1050 const internal::MatrixFreeFunctions::DoFInfo &
1059 template <
int dim,
typename Number>
1069 template <
int dim,
typename Number>
1081 template <
int dim,
typename Number>
1093 template <
int dim,
typename Number>
1095 std::pair<unsigned int,unsigned int>
1097 (
const std::pair<unsigned int,unsigned int> &range,
1098 const unsigned int degree,
1099 const unsigned int vector_component)
const 1102 if (
dof_info[vector_component].cell_active_fe_index.empty())
1105 if (
dof_info[vector_component].fe_index_conversion[0].first == degree)
1108 return std::pair<unsigned int,unsigned int> (range.second,range.second);
1111 const unsigned int fe_index =
1112 dof_info[vector_component].fe_index_from_degree(degree);
1113 if (fe_index >=
dof_info[vector_component].max_fe_index)
1114 return std::pair<unsigned int,unsigned int>(range.second, range.second);
1121 template <
int dim,
typename Number>
1123 std::pair<unsigned int,unsigned int>
1125 (
const std::pair<unsigned int,unsigned int> &range,
1126 const unsigned int fe_index,
1127 const unsigned int vector_component)
const 1130 const std::vector<unsigned int> &fe_indices =
1131 dof_info[vector_component].cell_active_fe_index;
1132 if (fe_indices.size() == 0)
1139 for (
unsigned int i=range.first+1; i<range.second; ++i)
1140 Assert (fe_indices[i] >= fe_indices[i-1],
1141 ExcMessage (
"Cell range must be over sorted range of fe indices in hp case!"));
1145 std::pair<unsigned int,unsigned int> return_range;
1146 return_range.first =
1147 std::lower_bound(&fe_indices[0] + range.first,
1148 &fe_indices[0] + range.second, fe_index)
1150 return_range.second =
1151 std::lower_bound(&fe_indices[0] + return_range.first,
1152 &fe_indices[0] + range.second,
1153 fe_index + 1)-&fe_indices[0];
1154 Assert(return_range.first >= range.first &&
1155 return_range.second <= range.second, ExcInternalError());
1156 return return_range;
1162 template <
int dim,
typename Number>
1166 const unsigned int vector_component)
1169 dof_info[vector_component].renumber_dofs (renumbering);
1174 template <
int dim,
typename Number>
1180 if (
dof_handlers.active_dof_handler == DoFHandlers::usual)
1188 Assert (
false, ExcNotImplemented());
1197 template <
int dim,
typename Number>
1201 const unsigned int vector_number,
1202 const unsigned int dof_index)
const 1209 const unsigned int irreg_filled =
dof_info[dof_index].row_starts[macro_cell_number][2];
1210 if (irreg_filled > 0)
1215 if (
dof_handlers.active_dof_handler == DoFHandlers::usual)
1224 "for underlying DoFHandler!"));
1227 std::pair<unsigned int,unsigned int> index =
1235 template <
int dim,
typename Number>
1239 const unsigned int vector_number,
1240 const unsigned int dof_index)
const 1247 const unsigned int irreg_filled =
dof_info[dof_index].row_starts[macro_cell_number][2];
1248 if (irreg_filled > 0)
1253 ExcNotImplemented());
1255 std::pair<unsigned int,unsigned int> index =
1263 template <
int dim,
typename Number>
1269 return dof_info[0].row_starts[macro_cell][2] > 0;
1274 template <
int dim,
typename Number>
1280 const unsigned int n_filled =
dof_info[0].row_starts[macro_cell][2];
1289 template <
int dim,
typename Number>
1293 const unsigned int active_fe_index)
const 1296 return dof_info[dof_index].dofs_per_cell[active_fe_index];
1301 template <
int dim,
typename Number>
1305 const unsigned int active_fe_index)
const 1309 return mapping_info.mapping_data_gen[quad_index].n_q_points[active_fe_index];
1314 template <
int dim,
typename Number>
1318 const unsigned int active_fe_index)
const 1321 return dof_info[dof_index].dofs_per_face[active_fe_index];
1326 template <
int dim,
typename Number>
1330 const unsigned int active_fe_index)
const 1334 return mapping_info.mapping_data_gen[quad_index].n_q_points_face[active_fe_index];
1339 template <
int dim,
typename Number>
1345 return dof_info[dof_index].vector_partitioner->locally_owned_range();
1350 template <
int dim,
typename Number>
1356 return dof_info[dof_index].vector_partitioner->ghost_indices();
1361 template <
int dim,
typename Number>
1365 const unsigned int index_quad,
1366 const unsigned int active_fe_index,
1367 const unsigned int active_quad_index)
const 1374 active_fe_index, active_quad_index);
1379 template <
int dim,
typename Number>
1383 const unsigned int active_fe_index)
const 1387 quadrature[active_fe_index];
1392 template <
int dim,
typename Number>
1396 const unsigned int active_fe_index)
const 1400 face_quadrature[active_fe_index];
1405 template <
int dim,
typename Number>
1415 template <
int dim,
typename Number>
1431 template <
typename DoFHandlerType>
1433 std::vector<IndexSet>
1434 extract_locally_owned_index_sets (
const std::vector<const DoFHandlerType *> &dofh,
1435 const unsigned int level)
1437 std::vector<IndexSet> locally_owned_set;
1438 locally_owned_set.reserve (dofh.size());
1439 for (
unsigned int j=0; j<dofh.size(); j++)
1441 locally_owned_set.push_back(dofh[j]->locally_owned_dofs());
1444 return locally_owned_set;
1447 template <
int dim,
int spacedim>
1449 std::vector<IndexSet>
1450 extract_locally_owned_index_sets (
const std::vector<const ::DoFHandler<dim,spacedim> *> &dofh,
1451 const unsigned int level)
1453 std::vector<IndexSet> locally_owned_set;
1454 locally_owned_set.reserve (dofh.size());
1455 for (
unsigned int j=0; j<dofh.size(); j++)
1457 locally_owned_set.push_back(dofh[j]->locally_owned_dofs());
1459 locally_owned_set.push_back(dofh[j]->locally_owned_mg_dofs(level));
1460 return locally_owned_set;
1467 template <
int dim,
typename Number>
1468 template <
typename DoFHandlerType,
typename QuadratureType>
1470 reinit(
const DoFHandlerType &dof_handler,
1472 const QuadratureType &quad,
1476 std::vector<const ConstraintMatrix *> constraints;
1477 std::vector<QuadratureType> quads;
1479 dof_handlers.push_back(&dof_handler);
1480 constraints.push_back (&constraints_in);
1481 quads.push_back (quad);
1483 std::vector<IndexSet> locally_owned_sets =
1484 internal::MatrixFree::extract_locally_owned_index_sets
1492 template <
int dim,
typename Number>
1493 template <
typename DoFHandlerType,
typename QuadratureType>
1496 const DoFHandlerType &dof_handler,
1498 const QuadratureType &quad,
1502 std::vector<const ConstraintMatrix *> constraints;
1503 std::vector<QuadratureType> quads;
1505 dof_handlers.push_back(&dof_handler);
1506 constraints.push_back (&constraints_in);
1507 quads.push_back (quad);
1509 std::vector<IndexSet> locally_owned_sets =
1510 internal::MatrixFree::extract_locally_owned_index_sets
1512 reinit(mapping, dof_handlers,constraints,locally_owned_sets, quads,
1518 template <
int dim,
typename Number>
1519 template <
typename DoFHandlerType,
typename QuadratureType>
1521 reinit(
const std::vector<const DoFHandlerType *> &dof_handler,
1522 const std::vector<const ConstraintMatrix *> &constraint,
1523 const std::vector<QuadratureType> &quad,
1526 std::vector<IndexSet> locally_owned_set =
1527 internal::MatrixFree::extract_locally_owned_index_sets
1530 static_cast<const std::vector<Quadrature<1>
>&>(quad),
1536 template <
int dim,
typename Number>
1537 template <
typename DoFHandlerType,
typename QuadratureType>
1539 reinit(
const std::vector<const DoFHandlerType *> &dof_handler,
1540 const std::vector<const ConstraintMatrix *> &constraint,
1541 const QuadratureType &quad,
1544 std::vector<QuadratureType> quads;
1545 quads.push_back(quad);
1546 std::vector<IndexSet> locally_owned_set =
1547 internal::MatrixFree::extract_locally_owned_index_sets
1555 template <
int dim,
typename Number>
1556 template <
typename DoFHandlerType,
typename QuadratureType>
1559 const std::vector<const DoFHandlerType *> &dof_handler,
1560 const std::vector<const ConstraintMatrix *> &constraint,
1561 const QuadratureType &quad,
1564 std::vector<QuadratureType> quads;
1565 quads.push_back(quad);
1566 std::vector<IndexSet> locally_owned_set =
1567 internal::MatrixFree::extract_locally_owned_index_sets
1569 reinit(mapping, dof_handler,constraint,locally_owned_set, quads,
1575 template <
int dim,
typename Number>
1576 template <
typename DoFHandlerType,
typename QuadratureType>
1579 const std::vector<const DoFHandlerType *> &dof_handler,
1580 const std::vector<const ConstraintMatrix *> &constraint,
1581 const std::vector<QuadratureType> &quad,
1584 std::vector<IndexSet> locally_owned_set =
1585 internal::MatrixFree::extract_locally_owned_index_sets
1587 reinit(mapping, dof_handler,constraint,locally_owned_set,
1588 quad, additional_data);
1593 template <
int dim,
typename Number>
1594 template <
typename DoFHandlerType,
typename QuadratureType>
1597 const std::vector<const DoFHandlerType *> &dof_handler,
1598 const std::vector<const ConstraintMatrix *> &constraint,
1599 const std::vector<IndexSet> &locally_owned_set,
1600 const std::vector<QuadratureType> &quad,
1604 std::vector<hp::QCollection<1> > quad_hp;
1605 for (
unsigned int q=0; q<quad.size(); ++q)
1609 constraint, locally_owned_set, quad_hp, additional_data);
1625 template<
typename VectorStruct>
1626 bool update_ghost_values_start_block (
const VectorStruct &vec,
1627 const unsigned int channel,
1629 template<
typename VectorStruct>
1630 void reset_ghost_values_block (
const VectorStruct &vec,
1631 const bool zero_out_ghosts,
1633 template<
typename VectorStruct>
1634 void update_ghost_values_finish_block (
const VectorStruct &vec,
1636 template<
typename VectorStruct>
1637 void compress_start_block (
const VectorStruct &vec,
1638 const unsigned int channel,
1640 template<
typename VectorStruct>
1641 void compress_finish_block (
const VectorStruct &vec,
1644 template<
typename VectorStruct>
1645 bool update_ghost_values_start_block (
const VectorStruct &,
1651 template<
typename VectorStruct>
1652 void reset_ghost_values_block (
const VectorStruct &,
1656 template<
typename VectorStruct>
1657 void update_ghost_values_finish_block (
const VectorStruct &,
1660 template<
typename VectorStruct>
1661 void compress_start_block (
const VectorStruct &,
1665 template<
typename VectorStruct>
1666 void compress_finish_block (
const VectorStruct &,
1674 template<
typename VectorStruct>
1676 bool update_ghost_values_start (
const VectorStruct &vec,
1677 const unsigned int channel = 0)
1680 update_ghost_values_start_block(vec, channel,
1686 template<
typename Number>
1689 const unsigned int channel = 0)
1693 return return_value;
1698 template <
typename VectorStruct>
1700 bool update_ghost_values_start (
const std::vector<VectorStruct> &vec)
1702 bool return_value =
false;
1703 for (
unsigned int comp=0; comp<vec.size(); comp++)
1704 return_value = update_ghost_values_start(vec[comp], comp);
1705 return return_value;
1710 template <
typename VectorStruct>
1712 bool update_ghost_values_start (
const std::vector<VectorStruct *> &vec)
1714 bool return_value =
false;
1715 for (
unsigned int comp=0; comp<vec.size(); comp++)
1716 return_value = update_ghost_values_start(*vec[comp], comp);
1717 return return_value;
1722 template<
typename VectorStruct>
1724 bool update_ghost_values_start_block (
const VectorStruct &vec,
1725 const unsigned int channel,
1728 bool return_value =
false;
1729 for (
unsigned int i=0; i<vec.n_blocks(); ++i)
1730 return_value = update_ghost_values_start(vec.block(i), channel+509*i);
1731 return return_value;
1739 template<
typename VectorStruct>
1741 void reset_ghost_values (
const VectorStruct &vec,
1742 const bool zero_out_ghosts)
1744 reset_ghost_values_block(vec, zero_out_ghosts,
1750 template<
typename Number>
1753 const bool zero_out_ghosts)
1755 if (zero_out_ghosts)
1761 template <
typename VectorStruct>
1763 void reset_ghost_values (
const std::vector<VectorStruct> &vec,
1764 const bool zero_out_ghosts)
1766 for (
unsigned int comp=0; comp<vec.size(); comp++)
1767 reset_ghost_values(vec[comp], zero_out_ghosts);
1772 template <
typename VectorStruct>
1774 void reset_ghost_values (
const std::vector<VectorStruct *> &vec,
1775 const bool zero_out_ghosts)
1777 for (
unsigned int comp=0; comp<vec.size(); comp++)
1778 reset_ghost_values(*vec[comp], zero_out_ghosts);
1783 template<
typename VectorStruct>
1785 void reset_ghost_values_block (
const VectorStruct &vec,
1786 const bool zero_out_ghosts,
1789 for (
unsigned int i=0; i<vec.n_blocks(); ++i)
1790 reset_ghost_values(vec.block(i), zero_out_ghosts);
1795 template <
typename VectorStruct>
1797 void update_ghost_values_finish (
const VectorStruct &vec)
1799 update_ghost_values_finish_block(vec,
1805 template <
typename Number>
1814 template <
typename VectorStruct>
1816 void update_ghost_values_finish (
const std::vector<VectorStruct> &vec)
1818 for (
unsigned int comp=0; comp<vec.size(); comp++)
1819 update_ghost_values_finish(vec[comp]);
1824 template <
typename VectorStruct>
1826 void update_ghost_values_finish (
const std::vector<VectorStruct *> &vec)
1828 for (
unsigned int comp=0; comp<vec.size(); comp++)
1829 update_ghost_values_finish(*vec[comp]);
1834 template <
typename VectorStruct>
1836 void update_ghost_values_finish_block (
const VectorStruct &vec,
1839 for (
unsigned int i=0; i<vec.n_blocks(); ++i)
1840 update_ghost_values_finish(vec.block(i));
1845 template <
typename VectorStruct>
1847 void compress_start (VectorStruct &vec,
1848 const unsigned int channel = 0)
1850 compress_start_block (vec, channel,
1856 template <
typename Number>
1859 const unsigned int channel = 0)
1866 template <
typename VectorStruct>
1868 void compress_start (std::vector<VectorStruct> &vec)
1870 for (
unsigned int comp=0; comp<vec.size(); comp++)
1871 compress_start (vec[comp], comp);
1876 template <
typename VectorStruct>
1878 void compress_start (std::vector<VectorStruct *> &vec)
1880 for (
unsigned int comp=0; comp<vec.size(); comp++)
1881 compress_start (*vec[comp], comp);
1886 template <
typename VectorStruct>
1888 void compress_start_block (VectorStruct &vec,
1889 const unsigned int channel,
1892 for (
unsigned int i=0; i<vec.n_blocks(); ++i)
1893 compress_start(vec.block(i), channel + 500*i);
1898 template <
typename VectorStruct>
1900 void compress_finish (VectorStruct &vec)
1902 compress_finish_block(vec,
1908 template <
typename Number>
1917 template <
typename VectorStruct>
1919 void compress_finish (std::vector<VectorStruct> &vec)
1921 for (
unsigned int comp=0; comp<vec.size(); comp++)
1922 compress_finish(vec[comp]);
1927 template <
typename VectorStruct>
1929 void compress_finish (std::vector<VectorStruct *> &vec)
1931 for (
unsigned int comp=0; comp<vec.size(); comp++)
1932 compress_finish(*vec[comp]);
1937 template <
typename VectorStruct>
1939 void compress_finish_block (VectorStruct &vec,
1942 for (
unsigned int i=0; i<vec.n_blocks(); ++i)
1943 compress_finish(vec.block(i));
1948 #ifdef DEAL_II_WITH_THREADS 1955 template<
typename Worker>
1956 class CellWork :
public tbb::task
1959 CellWork (
const Worker &worker_in,
1960 const unsigned int partition_in,
1962 const bool is_blocked_in)
1965 partition (partition_in),
1967 is_blocked (is_blocked_in)
1969 tbb::task *execute ()
1971 std::pair<unsigned int, unsigned int> cell_range
1972 (
task_info.partition_color_blocks_data[partition],
1973 task_info.partition_color_blocks_data[partition+1]);
1975 if (is_blocked==
true)
1976 dummy->spawn (*dummy);
1980 tbb::empty_task *dummy;
1983 const Worker &worker;
1984 const unsigned int partition;
1986 const bool is_blocked;
1991 template<
typename Worker>
1992 class PartitionWork :
public tbb::task
1995 PartitionWork (
const Worker &function_in,
1996 const unsigned int partition_in,
1998 const bool is_blocked_in =
false)
2000 function (function_in),
2001 partition (partition_in),
2003 is_blocked (is_blocked_in)
2005 tbb::task *execute ()
2007 tbb::empty_task *root =
new( tbb::task::allocate_root() )
2009 unsigned int evens =
task_info.partition_evens[partition];
2010 unsigned int odds =
task_info.partition_odds[partition];
2011 unsigned int n_blocked_workers =
2012 task_info.partition_n_blocked_workers[partition];
2013 unsigned int n_workers =
task_info.partition_n_workers[partition];
2014 std::vector<CellWork<Worker>*> worker(n_workers);
2015 std::vector<CellWork<Worker>*> blocked_worker(n_blocked_workers);
2017 root->set_ref_count(evens+1);
2018 for (
unsigned int j=0; j<evens; j++)
2020 worker[j] =
new(root->allocate_child())
2022 partition_color_blocks_row_index[partition]+2*j,
2026 worker[j]->set_ref_count(2);
2027 blocked_worker[j-1]->dummy =
new(worker[j]->allocate_child())
2029 worker[j-1]->spawn(*blocked_worker[j-1]);
2032 worker[j]->set_ref_count(1);
2035 blocked_worker[j] =
new(worker[j]->allocate_child())
2037 partition_color_blocks_row_index
2044 worker[evens] =
new(worker[j]->allocate_child())
2046 partition_color_blocks_row_index[partition]+2*j+1,
2048 worker[j]->spawn(*worker[evens]);
2052 tbb::empty_task *child =
new(worker[j]->allocate_child())
2054 worker[j]->spawn(*child);
2059 root->wait_for_all();
2060 root->destroy(*root);
2061 if (is_blocked==
true)
2062 dummy->spawn (*dummy);
2066 tbb::empty_task *dummy;
2069 const Worker &
function;
2070 const unsigned int partition;
2072 const bool is_blocked;
2081 template <
typename Worker>
2085 CellWork (
const Worker &worker_in,
2091 void operator()(
const tbb::blocked_range<unsigned int> &r)
const 2093 for (
unsigned int block=r.begin(); block<r.end(); block++)
2095 std::pair<unsigned int,unsigned int> cell_range;
2096 if (
task_info.position_short_block<block)
2098 cell_range.first = (block-1)*
task_info.block_size+
2100 cell_range.second = cell_range.first +
task_info.block_size;
2104 cell_range.first = block*
task_info.block_size;
2105 cell_range.second = cell_range.first +
2106 ((block ==
task_info.position_short_block)?
2109 worker (cell_range);
2113 const Worker &worker;
2118 template<
typename Worker>
2119 class PartitionWork :
public tbb::task
2122 PartitionWork (
const Worker &worker_in,
2123 const unsigned int partition_in,
2125 const bool is_blocked_in)
2128 partition (partition_in),
2130 is_blocked (is_blocked_in)
2132 tbb::task *execute ()
2134 unsigned int lower =
task_info.partition_color_blocks_data[partition],
2135 upper =
task_info.partition_color_blocks_data[partition+1];
2136 parallel_for(tbb::blocked_range<unsigned int>(lower,upper,1),
2138 if (is_blocked==
true)
2139 dummy->spawn (*dummy);
2143 tbb::empty_task *dummy;
2146 const Worker &worker;
2147 const unsigned int partition;
2149 const bool is_blocked;
2155 template<
typename VectorStruct>
2156 class MPIComDistribute :
public tbb::task
2159 MPIComDistribute (
const VectorStruct &src_in)
2164 tbb::task *execute ()
2166 internal::update_ghost_values_finish(src);
2171 const VectorStruct &src;
2176 template<
typename VectorStruct>
2177 class MPIComCompress :
public tbb::task
2180 MPIComCompress (VectorStruct &dst_in)
2185 tbb::task *execute ()
2187 internal::compress_start(dst);
2195 #endif // DEAL_II_WITH_THREADS 2201 template <
int dim,
typename Number>
2202 template <
typename OutVector,
typename InVector>
2209 const std::pair<
unsigned int,
2210 unsigned int> &)> &cell_operation,
2212 const InVector &src)
const 2215 bool ghosts_were_not_set = internal::update_ghost_values_start (src);
2217 #ifdef DEAL_II_WITH_THREADS 2226 std_cxx11::function<void (const std::pair<unsigned int,unsigned int> &range)>
2229 const Worker func = std_cxx11::bind (std_cxx11::ref(cell_operation),
2230 std_cxx11::cref(*
this),
2231 std_cxx11::ref(dst),
2232 std_cxx11::cref(src),
2235 if (
task_info.use_partition_partition ==
true)
2237 tbb::empty_task *root =
new( tbb::task::allocate_root() )
2241 root->set_ref_count(evens+1);
2242 unsigned int n_blocked_workers =
task_info.n_blocked_workers;
2243 unsigned int n_workers =
task_info.n_workers;
2244 std::vector<internal::partition::PartitionWork<Worker>*>
2246 std::vector<internal::partition::PartitionWork<Worker>*>
2247 blocked_worker(n_blocked_workers);
2248 internal::MPIComCompress<OutVector> *worker_compr =
2249 new(root->allocate_child())
2250 internal::MPIComCompress<OutVector>(dst);
2251 worker_compr->set_ref_count(1);
2252 for (
unsigned int j=0; j<evens; j++)
2256 worker[j] =
new(root->allocate_child())
2257 internal::partition::PartitionWork<Worker>
2259 worker[j]->set_ref_count(2);
2260 blocked_worker[j-1]->dummy =
new(worker[j]->allocate_child())
2263 worker[j-1]->spawn(*blocked_worker[j-1]);
2265 worker_compr->spawn(*blocked_worker[j-1]);
2269 worker[j] =
new(worker_compr->allocate_child())
2270 internal::partition::PartitionWork<Worker>
2272 worker[j]->set_ref_count(2);
2273 internal::MPIComDistribute<InVector> *worker_dist =
2274 new (worker[j]->allocate_child())
2275 internal::MPIComDistribute<InVector>(src);
2276 worker_dist->spawn(*worker_dist);
2280 blocked_worker[j] =
new(worker[j]->allocate_child())
2281 internal::partition::PartitionWork<Worker>
2288 worker[evens] =
new(worker[j]->allocate_child())
2289 internal::partition::PartitionWork<Worker>
2291 worker[j]->spawn(*worker[evens]);
2295 tbb::empty_task *child =
new(worker[j]->allocate_child())
2297 worker[j]->spawn(*child);
2302 root->wait_for_all();
2303 root->destroy(*root);
2314 tbb::empty_task *root =
new( tbb::task::allocate_root() ) tbb::empty_task;
2315 root->set_ref_count(evens+1);
2316 unsigned int n_blocked_workers = odds-(odds+evens+1)%2;
2317 unsigned int n_workers =
task_info.partition_color_blocks_data.size()-1-
2319 std::vector<internal::color::PartitionWork<Worker>*> worker(n_workers);
2320 std::vector<internal::color::PartitionWork<Worker>*> blocked_worker(n_blocked_workers);
2321 unsigned int worker_index = 0, slice_index = 0;
2322 unsigned int spawn_index = 0, spawn_index_new = 0;
2323 int spawn_index_child = -2;
2324 internal::MPIComCompress<OutVector> *worker_compr =
new(root->allocate_child())
2325 internal::MPIComCompress<OutVector>(dst);
2326 worker_compr->set_ref_count(1);
2327 for (
unsigned int part=0;
2328 part<
task_info.partition_color_blocks_row_index.size()-1; part++)
2330 spawn_index_new = worker_index;
2332 worker[worker_index] =
new(worker_compr->allocate_child())
2333 internal::color::PartitionWork<Worker>(func,slice_index,
task_info,
false);
2335 worker[worker_index] =
new(root->allocate_child())
2336 internal::color::PartitionWork<Worker>(func,slice_index,
task_info,
false);
2338 for (; slice_index<
task_info.partition_color_blocks_row_index[part+1];
2341 worker[worker_index]->set_ref_count(1);
2343 worker[worker_index] =
new (worker[worker_index-1]->allocate_child())
2344 internal::color::PartitionWork<Worker>(func,slice_index,
task_info,
false);
2346 worker[worker_index]->set_ref_count(2);
2349 blocked_worker[(part-1)/2]->dummy =
2350 new (worker[worker_index]->allocate_child()) tbb::empty_task;
2352 if (spawn_index_child == -1)
2353 worker[spawn_index]->spawn(*blocked_worker[(part-1)/2]);
2355 worker[spawn_index]->spawn(*worker[spawn_index_child]);
2356 spawn_index = spawn_index_new;
2357 spawn_index_child = -2;
2361 internal::MPIComDistribute<InVector> *worker_dist =
2362 new (worker[worker_index]->allocate_child())
2363 internal::MPIComDistribute<InVector>(src);
2364 worker_dist->spawn(*worker_dist);
2368 if (part<
task_info.partition_color_blocks_row_index.size()-1)
2370 if (part<
task_info.partition_color_blocks_row_index.size()-2)
2372 blocked_worker[part/2] =
new(worker[worker_index-1]->allocate_child())
2373 internal::color::PartitionWork<Worker>(func,slice_index,
task_info,
true);
2376 task_info.partition_color_blocks_row_index[part+1])
2378 blocked_worker[part/2]->set_ref_count(1);
2379 worker[worker_index] =
new(blocked_worker[part/2]->allocate_child())
2380 internal::color::PartitionWork<Worker>(func,slice_index,
task_info,
false);
2385 spawn_index_child = -1;
2389 for (; slice_index<
task_info.partition_color_blocks_row_index[part+1];
2393 task_info.partition_color_blocks_row_index[part])
2395 worker[worker_index]->set_ref_count(1);
2398 worker[worker_index] =
new (worker[worker_index-1]->allocate_child())
2399 internal::color::PartitionWork<Worker>(func,slice_index,
task_info,
false);
2401 spawn_index_child = worker_index;
2406 tbb::empty_task *
final =
new (worker[worker_index-1]->allocate_child())
2408 worker[spawn_index]->spawn(*
final);
2409 spawn_index_child = worker_index-1;
2413 worker[spawn_index]->spawn(*worker[spawn_index_child]);
2414 root->wait_for_all();
2415 root->destroy(*root);
2421 Assert(evens==1,ExcInternalError());
2422 internal::update_ghost_values_finish(src);
2424 for (
unsigned int color=0;
2425 color <
task_info.partition_color_blocks_row_index[1];
2428 unsigned int lower =
task_info.partition_color_blocks_data[color],
2429 upper =
task_info.partition_color_blocks_data[color+1];
2430 parallel_for(tbb::blocked_range<unsigned int>(lower,upper,1),
2431 internal::color::CellWork<Worker>
2435 internal::compress_start(dst);
2443 std::pair<unsigned int,unsigned int> cell_range;
2447 cell_range.first = 0;
2448 cell_range.second =
size_info.boundary_cells_start;
2449 cell_operation (*
this, dst, src, cell_range);
2454 internal::update_ghost_values_finish(src);
2459 cell_range.first =
size_info.boundary_cells_start;
2460 cell_range.second =
size_info.boundary_cells_end;
2461 cell_operation (*
this, dst, src, cell_range);
2464 internal::compress_start(dst);
2469 cell_range.first =
size_info.boundary_cells_end;
2470 cell_range.second =
size_info.n_macro_cells;
2471 cell_operation (*
this, dst, src, cell_range);
2476 internal::compress_finish(dst);
2477 internal::reset_ghost_values(src, ghosts_were_not_set);
2482 template <
int dim,
typename Number>
2483 template <
typename CLASS,
typename OutVector,
typename InVector>
2490 const std::pair<
unsigned int,
2491 unsigned int> &)
const,
2492 const CLASS *owning_class,
2494 const InVector &src)
const 2498 std_cxx11::function<void (const MatrixFree<dim,Number> &,
2501 const std::pair<
unsigned int,
2503 function = std_cxx11::bind<void>(function_pointer,
2514 template <
int dim,
typename Number>
2515 template <
typename CLASS,
typename OutVector,
typename InVector>
2522 const std::pair<
unsigned int,
2524 CLASS *owning_class,
2526 const InVector &src)
const 2530 std_cxx11::function<void (const MatrixFree<dim,Number> &,
2533 const std::pair<
unsigned int,
2535 function = std_cxx11::bind<void>(function_pointer,
2545 #endif // ifndef DOXYGEN 2549 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info(const unsigned int fe_component=0, const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0, const unsigned int hp_active_quad_index=0) const
const Quadrature< dim-1 > & get_face_quadrature(const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
unsigned int n_physical_cells() const
DoFHandler< dim >::cell_iterator get_cell_iterator(const unsigned int macro_cell_number, const unsigned int vector_number, const unsigned int fe_component=0) const
void reinit(const Mapping< dim > &mapping, const DoFHandlerType &dof_handler, const ConstraintMatrix &constraint, const IndexSet &locally_owned_dofs, const QuadratureType &quad, const AdditionalData additional_data=AdditionalData())
internal::MatrixFreeFunctions::TaskInfo task_info
bool indices_are_initialized
static const unsigned int invalid_unsigned_int
unsigned int get_dofs_per_face(const unsigned int fe_component=0, const unsigned int hp_active_fe_index=0) const
const Number * constraint_pool_begin(const unsigned int pool_index) const
std::vector< unsigned int > constraint_pool_row_index
const DoFHandler< dim > & get_dof_handler(const unsigned int fe_component=0) const
unsigned int get_n_q_points_face(const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
#define AssertDimension(dim1, dim2)
AdditionalData(const MPI_Comm mpi_communicator=MPI_COMM_SELF, const TasksParallelScheme tasks_parallel_scheme=partition_partition, const unsigned int tasks_block_size=0, const UpdateFlags mapping_update_flags=update_gradients|update_JxW_values, const unsigned int level_mg_handler=numbers::invalid_unsigned_int, const bool store_plain_indices=true, const bool initialize_indices=true, const bool initialize_mapping=true)
const internal::MatrixFreeFunctions::DoFInfo & get_dof_info(const unsigned int fe_component=0) const
internal::MatrixFreeFunctions::MappingInfo< dim, Number > mapping_info
unsigned int n_components_filled(const unsigned int macro_cell_number) const
unsigned int n_components() const
bool has_ghost_elements() const
::ExceptionBase & ExcMessage(std::string arg1)
std::vector< Number > constraint_pool_data
#define AssertIndexRange(index, range)
Table< 4, internal::MatrixFreeFunctions::ShapeInfo< Number > > shape_info
void copy_from(const MatrixFree< dim, Number > &matrix_free_base)
std::vector< std::pair< unsigned int, unsigned int > > cell_level_index
#define AssertThrow(cond, exc)
const std::vector< unsigned int > & get_constrained_dofs(const unsigned int fe_component=0) const
void compress_finish(::VectorOperation::values operation)
const IndexSet & get_ghost_set(const unsigned int fe_component=0) const
std::pair< unsigned int, unsigned int > create_cell_subrange_hp(const std::pair< unsigned int, unsigned int > &range, const unsigned int fe_degree, const unsigned int vector_component=0) const
ActiveSelector::active_cell_iterator active_cell_iterator
bool mapping_is_initialized
internal::MatrixFreeFunctions::SizeInfo size_info
ActiveSelector::cell_iterator cell_iterator
bool at_irregular_cell(const unsigned int macro_cell_number) const
#define Assert(cond, exc)
void print_memory_consumption(StreamType &out) const
unsigned int get_dofs_per_cell(const unsigned int fe_component=0, const unsigned int hp_active_fe_index=0) const
const Number * constraint_pool_end(const unsigned int pool_index) const
unsigned int n_constraint_pool_entries() const
unsigned int n_macro_cells() const
void initialize_dof_handlers(const std::vector< const DoFHandler< dim > *> &dof_handlers, const unsigned int level)
ActiveSelector::cell_iterator cell_iterator
void internal_reinit(const Mapping< dim > &mapping, const std::vector< const DoFHandler< dim > *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const std::vector< IndexSet > &locally_owned_set, const std::vector< hp::QCollection< 1 > > &quad, const AdditionalData additional_data)
TasksParallelScheme tasks_parallel_scheme
const IndexSet & get_locally_owned_set(const unsigned int fe_component=0) const
bool mapping_initialized() const
MPI_Comm mpi_communicator
std::pair< unsigned int, unsigned int > create_cell_subrange_hp_by_index(const std::pair< unsigned int, unsigned int > &range, const unsigned int fe_index, const unsigned int vector_component=0) const
void update_ghost_values_finish() const
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void cell_loop(const std_cxx11::function< void(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, OutVector &dst, const InVector &src) const
const internal::MatrixFreeFunctions::SizeInfo & get_size_info() const
unsigned int tasks_block_size
const std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > & get_vector_partitioner(const unsigned int vector_component=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
unsigned int level_mg_handler
hp::DoFHandler< dim >::active_cell_iterator get_hp_cell_iterator(const unsigned int macro_cell_number, const unsigned int vector_number, const unsigned int fe_component=0) const
Shape function gradients.
unsigned int get_n_q_points(const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
UpdateFlags mapping_update_flags
std::size_t memory_consumption() const
const Triangulation< dim, spacedim > & get_triangulation() const
void initialize_dof_vector(VectorType &vec, const unsigned int vector_component=0) const
bool indices_initialized() const
unsigned int size(const unsigned int i) const
void compress_start(const unsigned int communication_channel=0, ::VectorOperation::values operation=VectorOperation::add)
const Quadrature< dim > & get_quadrature(const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
std::vector< internal::MatrixFreeFunctions::DoFInfo > dof_info
void update_ghost_values_start(const unsigned int communication_channel=0) const
void renumber_dofs(std::vector< types::global_dof_index > &renumbering, const unsigned int vector_component=0)
const internal::MatrixFreeFunctions::TaskInfo & get_task_info() const
void print(std::ostream &out) const