Epetra Package Browser (Single Doxygen Collection)  Development
example/UG_Ex1/cxx_main.cpp
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ************************************************************************
4 //
5 // Epetra: Linear Algebra Services Package
6 // Copyright 2011 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ************************************************************************
41 //@HEADER
42 */
43 
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include "Epetra_config.h"
47 #ifdef HAVE_MPI
48 #include "mpi.h"
49 #include "Epetra_MpiComm.h"
50 #else
51 #include "Epetra_SerialComm.h"
52 #endif
53 #include "Epetra_Map.h"
56 #include "Epetra_Vector.h"
57 #include "Epetra_CrsMatrix.h"
58 #include "Epetra_Version.h"
59 
60 // prototype
61 double power_method(const Epetra_CrsMatrix& A);
62 int main(int argc, char *argv[]) {
63 #ifdef HAVE_MPI
64  MPI_Init(&argc,&argv);
65  Epetra_MpiComm Comm(MPI_COMM_WORLD);
66 #else
67  Epetra_SerialComm Comm;
68 #endif
69 
70  if (Comm.MyPID()==0)
71  std::cout << Epetra_Version() << std::endl << std::endl;
72 
73  std::cout << Comm << std::endl; // Print out process information
74  // Get the number of global equations from the command line
75  if (argc!=2) {
76  std::cout << "Usage: " << argv[0] << " number_of_equations" << std::endl;
77  exit(1);
78  }
79  int NumGlobalElements = atoi(argv[1]);
80  // Construct a Map that puts approximately the same number of
81  // equations on each processor.
82  Epetra_Map Map(NumGlobalElements, 0, Comm);
83  // Get update list and number of local equations from newly created Map.
84  int NumMyElements = Map.NumMyElements();
85  // Create an Epetra_CrsMatrix
86  Epetra_CrsMatrix A(Copy, Map, 1);
87  // Add rows one-at-a-time
88  Epetra_SerialDenseVector DiagVal(1);
89  DiagVal[0] = 2.0; // We set all diagonal values to 2
90  Epetra_IntSerialDenseVector ColIndices(1);
91  for (int i=0; i<NumMyElements; i++) {
92  int RowIndex = Map.GID(i);
93  ColIndices[0] = RowIndex;
94  // Put in the diagonal entry
95  A.InsertGlobalValues(RowIndex, DiagVal.Length(),
96  DiagVal.Values(), ColIndices.Values());
97  }
98  if (Map.MyGID(0)) { // Change the first global diagonal value to 4.0
99  DiagVal[0] = 4.0;
100  int RowIndex = 0;
101  ColIndices[0] = RowIndex;
102  A.ReplaceGlobalValues(RowIndex, DiagVal.Length(),
103  DiagVal.Values(), ColIndices.Values());
104  }
105  // Finish up
106  A.FillComplete();
107  // Iterate
108  double lambda = power_method(A);
109  if (Comm.MyPID()==0)
110  std::cout << std::endl << "Estimate of Dominant Eigenvalue = " << lambda << std::endl;
111 #ifdef UG_EX1_MPI
112  MPI_Finalize() ;
113 #endif
114 return 0;
115 }
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
Epetra_IntSerialDenseVector: A class for constructing and using dense vectors.
double * Values() const
Returns pointer to the values in vector.
Epetra_MpiComm: The Epetra MPI Communication Class.
std::string Epetra_Version()
Epetra_SerialDenseVector: A class for constructing and using dense vectors.
bool MyGID(int GID_in) const
Returns true if the GID passed in belongs to the calling processor in this map, otherwise returns fal...
int NumMyElements() const
Number of elements on the calling processor.
double power_method(const Epetra_CrsMatrix &A)
Epetra_SerialComm: The Epetra Serial Communication Class.
int MyPID() const
Return my process ID.
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed r...
int main(int argc, char *argv[])
int Length() const
Returns length of vector.
int * Values()
Returns pointer to the values in vector.
int GID(int LID) const
Returns global ID of local ID, return IndexBase-1 if not found on this processor. ...