Galeri Package Browser (Single Doxygen Collection)
Version of the Day
Main Page
Related Pages
Namespaces
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
example
LinearProblem.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_Maps.h"
11
#include "Galeri_CrsMatrices.h"
12
#include "Galeri_Utils.h"
13
#ifdef HAVE_MPI
14
#include "
Epetra_MpiComm.h
"
15
#include "mpi.h"
16
#else
17
#include "
Epetra_SerialComm.h
"
18
#endif
19
#include "
Epetra_Map.h
"
20
#include "
Epetra_Vector.h
"
21
#include "
Epetra_CrsMatrix.h
"
22
#include "
Epetra_LinearProblem.h
"
23
#include "
Teuchos_ParameterList.hpp
"
24
25
using namespace
Galeri;
26
27
// =========== //
28
// main driver //
29
// =========== //
30
31
int
main
(
int
argc,
char
* argv[])
32
{
33
#ifdef HAVE_MPI
34
MPI_Init(&argc, &argv);
35
Epetra_MpiComm
Comm(MPI_COMM_WORLD);
36
#else
37
Epetra_SerialComm
Comm;
38
#endif
39
40
// Here we create the linear problem
41
//
42
// Matrix * LHS = RHS
43
//
44
// with Matrix arising from a 5-point formula discretization.
45
46
Epetra_Map
* Map = 0;
47
Epetra_RowMatrix
* Matrix = 0;
48
49
Teuchos::ParameterList
GaleriList;
50
// dimension of the problem is nx x ny
51
GaleriList.
set
(
"nx"
, 10 * Comm.
NumProc
());
52
GaleriList.
set
(
"ny"
, 10);
53
// total number of processors is mx x my
54
GaleriList.
set
(
"mx"
, Comm.
NumProc
());
55
GaleriList.
set
(
"my"
, 1);
56
57
try
58
{
59
#ifndef GALERI_TEST_USE_LONGLONG_GO
60
Map = CreateMap(
"Cartesian2D"
, Comm, GaleriList);
61
#else
62
Map = CreateMap64(
"Cartesian2D"
, Comm, GaleriList);
63
#endif
64
Matrix = CreateCrsMatrix(
"Laplace2D"
, Map, GaleriList);
65
Epetra_Vector
ExactSolution(*Map); ExactSolution.Random();
66
Epetra_Vector
LHS(*Map); LHS.PutScalar(0.0);
67
Epetra_Vector
RHS(*Map);
68
69
Matrix->
Multiply
(
false
, ExactSolution, RHS);
70
71
Epetra_LinearProblem
Problem(Matrix, &LHS, &RHS);
72
73
// at this point any object that understand Epetra_LinearProblem can be
74
// used, for example AztecOO, Amesos. IFPACK and ML can be used to define a
75
// preconditioner for Matrix. Here we use a simple solver, based on
76
// LAPACK, that is meant for simple testing only.
77
78
Solve
(Problem);
79
80
// and we compute the norm of the true residual.
81
double
ResidualNorm = ComputeNorm(Matrix, &LHS, &RHS);
82
83
if
(Comm.
MyPID
() == 0)
84
cout << ResidualNorm << endl;
85
86
delete
Map;
87
delete
Matrix;
88
}
89
catch
(Galeri::Exception& rhs)
90
{
91
if
(Comm.
MyPID
() == 0)
92
{
93
cerr <<
"Caught exception: "
;
94
rhs.Print();
95
}
96
}
97
98
#ifdef HAVE_MPI
99
MPI_Finalize();
100
#endif
101
102
return
(EXIT_SUCCESS);
103
}
Epetra_Map
Epetra_LinearProblem.h
Epetra_Vector
Teuchos::ParameterList::set
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Epetra_SerialComm::MyPID
int MyPID() const
Epetra_MpiComm
Epetra_SerialComm.h
Teuchos_ParameterList.hpp
Epetra_SerialComm::NumProc
int NumProc() const
Teuchos::ParameterList
Epetra_SerialComm
Epetra_MpiComm.h
Epetra_RowMatrix::Multiply
virtual int Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
main
int main(int argc, char *argv[])
Definition:
CrsMatrix.cpp:29
Epetra_Map.h
Epetra_Vector.h
Solve
int Solve(int, TYPE *, TYPE *, TYPE *)
Epetra_LinearProblem
Epetra_RowMatrix
Epetra_CrsMatrix.h
Generated on Fri Dec 27 2024 09:26:21 for Galeri Package Browser (Single Doxygen Collection) by
1.8.5