43 #ifndef IFPACK2_CREATEOVERLAPGRAPH_HPP 
   44 #define IFPACK2_CREATEOVERLAPGRAPH_HPP 
   46 #include "Ifpack2_ConfigDefs.hpp" 
   47 #include "Tpetra_CrsGraph.hpp" 
   48 #include "Tpetra_CrsMatrix.hpp" 
   49 #include "Tpetra_Import.hpp" 
   50 #include "Teuchos_RefCountPtr.hpp" 
   71 template<
class GraphType>
 
   74                     const int overlapLevel)
 
   78   typedef Tpetra::Map<
typename GraphType::local_ordinal_type,
 
   79                       typename GraphType::global_ordinal_type,
 
   80                       typename GraphType::node_type> map_type;
 
   81   typedef Tpetra::Import<
typename GraphType::local_ordinal_type,
 
   82                          typename GraphType::global_ordinal_type,
 
   83                          typename GraphType::node_type> import_type;
 
   85     overlapLevel < 0, std::invalid_argument,
 
   86     "Ifpack2::createOverlapGraph: overlapLevel must be >= 0, " 
   87     "but you specified overlapLevel = " << overlapLevel << 
".");
 
   89   const int numProcs = inputGraph->getMap ()->getComm ()->getSize ();
 
   90   if (overlapLevel == 0 || numProcs < 2) {
 
   94   RCP<const map_type> overlapRowMap = inputGraph->getRowMap ();
 
   95   RCP<const map_type> domainMap = inputGraph->getDomainMap ();
 
   96   RCP<const map_type> rangeMap = inputGraph->getRangeMap ();
 
   98   RCP<GraphType> overlapGraph;
 
   99   RCP<const GraphType> oldGraph;
 
  100   RCP<const map_type> oldRowMap;
 
  101   for (
int level = 0; level < overlapLevel; ++level) {
 
  102     oldGraph = overlapGraph;
 
  103     oldRowMap = overlapRowMap;
 
  105     RCP<const import_type> overlapImporter;
 
  107       overlapImporter = inputGraph->getImporter ();
 
  109       overlapImporter = oldGraph->getImporter ();
 
  112     overlapRowMap = overlapImporter->getTargetMap ();
 
  113     if (level < overlapLevel - 1) {
 
  114       overlapGraph = 
rcp (
new GraphType (overlapRowMap, 0));
 
  119       overlapGraph = 
rcp (
new GraphType (overlapRowMap, overlapRowMap, 0));
 
  122     overlapGraph->doImport (*inputGraph, *overlapImporter, Tpetra::INSERT);
 
  123     overlapGraph->fillComplete (domainMap, rangeMap);
 
  138 template<
class MatrixType>
 
  141                      const int overlapLevel)
 
  145   typedef typename MatrixType::map_type map_type;
 
  146   typedef Tpetra::Import<
typename MatrixType::local_ordinal_type,
 
  147     typename MatrixType::global_ordinal_type,
 
  148     typename MatrixType::node_type> import_type;
 
  151     overlapLevel < 0, std::invalid_argument,
 
  152     "Ifpack2::createOverlapMatrix: overlapLevel must be >= 0, " 
  153     "but you specified overlapLevel = " << overlapLevel << 
".");
 
  155   const int numProcs = inputMatrix->getMap ()->getComm ()->getSize ();
 
  156   if (overlapLevel == 0 || numProcs < 2) {
 
  160   RCP<const map_type> overlapRowMap = inputMatrix->getRowMap ();
 
  161   RCP<const map_type> domainMap = inputMatrix->getDomainMap ();
 
  162   RCP<const map_type> rangeMap = inputMatrix->getRangeMap ();
 
  164   RCP<MatrixType> overlapMatrix;
 
  165   RCP<const MatrixType> oldMatrix;
 
  166   RCP<const map_type> oldRowMap;
 
  167   for (
int level = 0; level < overlapLevel; ++level) {
 
  168     oldMatrix = overlapMatrix;
 
  169     oldRowMap = overlapRowMap;
 
  171     RCP<const import_type> overlapImporter;
 
  173       overlapImporter = inputMatrix->getGraph ()->getImporter ();
 
  175       overlapImporter = oldMatrix->getGraph ()->getImporter ();
 
  178     overlapRowMap = overlapImporter->getTargetMap ();
 
  179     if (level < overlapLevel - 1) {
 
  180       overlapMatrix = 
rcp (
new MatrixType (overlapRowMap, 0));
 
  185       overlapMatrix = 
rcp (
new MatrixType (overlapRowMap, overlapRowMap, 0));
 
  188     overlapMatrix->doImport (*inputMatrix, *overlapImporter, Tpetra::INSERT);
 
  189     overlapMatrix->fillComplete (domainMap, rangeMap);
 
  192   return overlapMatrix;
 
  197 #endif // IFPACK2_CREATEOVERLAPGRAPH_HPP 
Teuchos::RCP< const GraphType > createOverlapGraph(const Teuchos::RCP< const GraphType > &inputGraph, const int overlapLevel)
Construct an overlapped graph for use with Ifpack2 preconditioners. 
Definition: Ifpack2_CreateOverlapGraph.hpp:73
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
Teuchos::RCP< const MatrixType > createOverlapMatrix(const Teuchos::RCP< const MatrixType > &inputMatrix, const int overlapLevel)
Construct an overlapped matrix for use with Ifpack2 preconditioners. 
Definition: Ifpack2_CreateOverlapGraph.hpp:140