xref: /petsc/src/sys/tests/ex45.cxx (revision 1bb3edfd3d106a9ce6e99783e45c9200c8e18d90)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Demonstrates call PETSc first and then Trilinos in the same program.\n\n";
3c4762a1bSJed Brown 
4*1bb3edfdSPatrick Sanan /*
5c4762a1bSJed Brown    Example obtained from: http://trilinos.org/docs/dev/packages/tpetra/doc/html/Tpetra_Lesson01.html
6*1bb3edfdSPatrick Sanan */
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscsys.h>
9c4762a1bSJed Brown #include <Teuchos_DefaultMpiComm.hpp> // wrapper for MPI_Comm
10c4762a1bSJed Brown #include <Tpetra_Version.hpp> // Tpetra version string
11c4762a1bSJed Brown 
12c4762a1bSJed Brown // Do something with the given communicator.  In this case, we just
13c4762a1bSJed Brown // print Tpetra's version to stdout on Process 0 in the given
14c4762a1bSJed Brown // communicator.
15c4762a1bSJed Brown void
16c4762a1bSJed Brown exampleRoutine (const Teuchos::RCP<const Teuchos::Comm<int> >& comm)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   if (comm->getRank () == 0) {
19c4762a1bSJed Brown     // On (MPI) Process 0, print out the Tpetra software version.
20c4762a1bSJed Brown     std::cout << Tpetra::version () << std::endl << std::endl;
21c4762a1bSJed Brown   }
22c4762a1bSJed Brown }
23c4762a1bSJed Brown 
24c4762a1bSJed Brown int main(int argc,char **argv)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   // These "using" declarations make the code more concise, in that
27c4762a1bSJed Brown   // you don't have to write the namespace along with the class or
28c4762a1bSJed Brown   // object name.  This is especially helpful with commonly used
29c4762a1bSJed Brown   // things like std::endl or Teuchos::RCP.
30c4762a1bSJed Brown   using std::cout;
31c4762a1bSJed Brown   using std::endl;
32c4762a1bSJed Brown   using Teuchos::Comm;
33c4762a1bSJed Brown   using Teuchos::MpiComm;
34c4762a1bSJed Brown   using Teuchos::RCP;
35c4762a1bSJed Brown   using Teuchos::rcp;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown     Every PETSc routine should begin with the PetscInitialize() routine.
39c4762a1bSJed Brown     argc, argv - These command line arguments are taken to extract the options
40c4762a1bSJed Brown                  supplied to PETSc and options supplied to MPI.
41c4762a1bSJed Brown     help       - When PETSc executable is invoked with the option -help,
42c4762a1bSJed Brown                  it prints the various options that can be applied at
43c4762a1bSJed Brown                  runtime.  The user can use the "help" variable place
44c4762a1bSJed Brown                  additional help messages in this printout.
45c4762a1bSJed Brown   */
46b8abcfdeSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
47c4762a1bSJed Brown   RCP<const Comm<int> > comm (new MpiComm<int> (PETSC_COMM_WORLD));
48c4762a1bSJed Brown   // Get my process' rank, and the total number of processes.
49c4762a1bSJed Brown   // Equivalent to MPI_Comm_rank resp. MPI_Comm_size.
50c4762a1bSJed Brown   const int myRank = comm->getRank ();
51c4762a1bSJed Brown   const int size = comm->getSize ();
52c4762a1bSJed Brown   if (myRank == 0) {
53c4762a1bSJed Brown     cout << "Total number of processes: " << size << endl;
54c4762a1bSJed Brown   }
55c4762a1bSJed Brown   // Do something with the new communicator.
56c4762a1bSJed Brown   exampleRoutine (comm);
57c4762a1bSJed Brown   // This tells the Trilinos test framework that the test passed.
58c4762a1bSJed Brown   if (myRank == 0) {
59c4762a1bSJed Brown     cout << "End Result: TEST PASSED" << endl;
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown 
62b8abcfdeSJacob Faibussowitsch   PetscCall(PetscFinalize());
63b8abcfdeSJacob Faibussowitsch   return 0;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /*TEST
67c4762a1bSJed Brown 
68c4762a1bSJed Brown    build:
69c4762a1bSJed Brown      requires: trilinos
70c4762a1bSJed Brown 
71c4762a1bSJed Brown    test:
72c4762a1bSJed Brown       nsize: 3
73c4762a1bSJed Brown       filter: grep -v "Tpetra in Trilinos"
74c4762a1bSJed Brown 
75c4762a1bSJed Brown TEST*/
76