Galeri Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
cxx_main_tpetra.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Galeri: Finite Element and Matrix Generation Package
4 //
5 // Copyright 2006 ETHZ/NTESS and the Galeri contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #include "Galeri_XpetraMaps.hpp"
11 #include "Galeri_MatrixTraits.hpp"
12 #include "Galeri_XpetraMatrixTypes.hpp"
13 #include "Galeri_XpetraProblemFactory.hpp"
14 
15 #include "Teuchos_DefaultComm.hpp"
17 
18 #ifndef GALERI_TEST_USE_LONGLONG_GO
19 #define GO int
20 #else
21 #define GO long long
22 #endif
23 #define Scalar int
24 #define LO int
25 #define Node Tpetra::KokkosClassic::DefaultNode::DefaultNodeType
26 
27 using namespace Galeri;
28 
29 int main(int argc, char* argv[])
30 {
31  using Teuchos::RCP;
32  using Teuchos::rcp;
33 
34  typedef Tpetra::Map<LO, GO, Node> Tpetra_Map;
35  typedef Tpetra::CrsMatrix<Scalar, LO, GO, Node> Tpetra_CrsMatrix;
36  typedef Tpetra::MultiVector<Scalar, LO, GO, Node> Tpetra_MultiVector;
37  typedef Teuchos::ScalarTraits<Scalar> ScalarTraits;
38 
39 #ifdef HAVE_MPI
40  MPI_Init(&argc, &argv);
41 #endif
42 
43  // Create comm
44  RCP<const Teuchos::Comm<int>> comm = Teuchos::DefaultComm<int>::getComm();
45 
46  // Here we create the linear problem
47  //
48  // Matrix * LHS = RHS
49  //
50  // with Matrix arising from a 5-point formula discretization.
51 
52  std::string mapType = "Cartesian2D";
53  auto mapParameters = Teuchos::ParameterList("Tpetra::Map");
54  // dimension of the problem is nx x ny
55  mapParameters.set("nx", 10 * comm->getSize());
56  mapParameters.set("ny", 10);
57  // total number of processors is mx x my
58  mapParameters.set("mx", comm->getSize());
59  mapParameters.set("my", 1);
60 
61  auto out = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
62 
63  try
64  {
65  // Creation of the map
66  auto map = RCP{Galeri::Xpetra::CreateMap<Scalar, GO, Tpetra_Map>(mapType, comm, mapParameters)};
67 
68  // Creation of linear problem
69  auto problem = Galeri::Xpetra::BuildProblem<Scalar, LO, GO, Tpetra_Map, Tpetra_CrsMatrix, Tpetra_MultiVector>("Laplace2D", map, mapParameters);
70 
71  // Build Matrix and MultiVectors
72  auto matrix = problem->BuildMatrix();
73  auto LHS = rcp(new Tpetra_MultiVector(matrix->getDomainMap(), 1));
74  auto RHS = rcp(new Tpetra_MultiVector(matrix->getRangeMap(), 1));
75  auto ExactSolution = rcp(new Tpetra_MultiVector(matrix->getDomainMap(), 1));
76 
77  ExactSolution->randomize(0, 100);
78  LHS->putScalar(ScalarTraits::zero());
79 
80  matrix->apply(*ExactSolution, *RHS);
81 
82  matrix->describe(*out, Teuchos::EVerbosityLevel::VERB_EXTREME);
83  LHS->describe(*out, Teuchos::EVerbosityLevel::VERB_EXTREME);
84  RHS->describe(*out, Teuchos::EVerbosityLevel::VERB_EXTREME);
85  ExactSolution->describe(*out, Teuchos::EVerbosityLevel::VERB_EXTREME);
86 
87  // at this point any LinearSolver can be used which understands the Tpetra objects. For example: Amesos2 or Ifpack2
88  }
89  catch (Galeri::Exception &rhs)
90  {
91  if (comm->getRank() == 0)
92  {
93  cerr << "Caught exception: ";
94  rhs.Print();
95 
96 #ifdef HAVE_MPI
97  MPI_Finalize();
98 #endif
99  return (EXIT_FAILURE);
100  }
101  }
102 
103 #ifdef HAVE_MPI
104  MPI_Finalize();
105 #endif
106 
107  return (EXIT_SUCCESS);
108 }
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int main(int argc, char *argv[])
Definition: CrsMatrix.cpp:29