LifeV
lifev/core/testsuite/vector_container/main.cpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28  * @file
29  * @brief File containing the vector container test
30  *
31  * @date 29-09-2009
32  * @author Cristiano Malossi <cristiano.malossi@epfl.ch>
33  *
34  * @maintainer Cristiano Malossi <cristiano.malossi@epfl.ch>
35  *
36  * This is a test to verify that the VectorContainer works correctly.
37  */
38 
39 
40 
41 #include <Epetra_ConfigDefs.h>
42 #ifdef HAVE_MPI
43 #include <mpi.h>
44 #include <Epetra_MpiComm.h>
45 #else
46 #include <Epetra_SerialComm.h>
47 #endif
48 
49 
50 #include <lifev/core/array/MapEpetra.hpp>
51 #include <lifev/core/array/VectorEpetra.hpp>
52 
53 #include "TestFunction.hpp"
54 
55 using namespace LifeV;
56 
57 Int
58 main ( Int argc, char** argv )
59 {
60 #ifdef HAVE_MPI
61  std::cout << "MPI Initialization" << std::endl;
62  MPI_Init ( &argc, &argv );
63 #endif
64 
65  // ===================================================
66  // TEST: LIFEV::EPETRAVECTOR
67  // ===================================================
68 
69  //Setup main communicator
70  std::shared_ptr<Epetra_Comm> comm;
71 
72  Int nprocs;
73  Int rank;
74 
75  //MPI Preprocessing
76 #ifdef EPETRA_MPI
77 
78  MPI_Comm_size ( MPI_COMM_WORLD, &nprocs );
79  MPI_Comm_rank ( MPI_COMM_WORLD, &rank );
80 
81  if ( rank == 0 )
82  {
83  std::cout << "MPI processes: " << nprocs << std::endl;
84 
85  std::cout << "MPI Epetra Initialization ... " << std::endl;
86  }
87  comm.reset ( new Epetra_MpiComm ( MPI_COMM_WORLD ) );
88 
89  comm->Barrier();
90 
91 #else
92 
93  std::cout << "MPI SERIAL Epetra Initialization ... " << std::endl;
94  comm.reset ( new Epetra_SerialComm() );
95 #endif
96 
97  std::cout << std::setprecision (15) << std::endl;
98 
99  //BASE VECTORS
100  typedef VectorEpetra Vector;
101  typedef std::shared_ptr<Vector> Vector_ptr;
102 
103  Int MyGlobalIElementsA[3], MyGlobalIElementsB[2];
104  MyGlobalIElementsA[0] = 0;
105  MyGlobalIElementsA[1] = 1;
106  MyGlobalIElementsA[2] = 2;
107  MyGlobalIElementsB[0] = 0;
108  MyGlobalIElementsB[1] = 1;
109 
110  MapEpetra mapA ( 3, 3, &MyGlobalIElementsA[0], comm);
111  MapEpetra mapB ( 2, 2, &MyGlobalIElementsB[0], comm);
112 
113  Vector_ptr A1, B1, A2, B2, A3, B3, A4, B4;
114  A1.reset ( new Vector ( mapA, Unique ) );
115  B1.reset ( new Vector ( mapB, Unique ) );
116  A2.reset ( new Vector ( mapA, Unique ) );
117  B2.reset ( new Vector ( mapB, Unique ) );
118  A3.reset ( new Vector ( mapA, Unique ) );
119  B3.reset ( new Vector ( mapB, Unique ) );
120  A4.reset ( new Vector ( mapA, Unique ) );
121  B4.reset ( new Vector ( mapB, Unique ) );
122 
123  //Fill the vectors
124  for ( UInt i (0) ; i < static_cast<UInt> (A1->size() ) ; ++i )
125  {
126  (*A1) [i] = i;
127  }
128  for ( UInt i (0) ; i < static_cast<UInt> (B1->size() ) ; ++i )
129  {
130  (*B1) [i] = A1->size() + i;
131  }
132 
133  for ( UInt i (0) ; i < static_cast<UInt> (A2->size() ) ; ++i )
134  {
135  (*A2) [i] = i * 10;
136  }
137  for ( UInt i (0) ; i < static_cast<UInt> (B2->size() ) ; ++i )
138  {
139  (*B2) [i] = (A2->size() + i) * 10;
140  }
141 
142  for ( UInt i (0) ; i < static_cast<UInt> (A3->size() ) ; ++i )
143  {
144  (*A3) [i] = i * 10;
145  }
146  for ( UInt i (0) ; i < static_cast<UInt> (B3->size() ) ; ++i )
147  {
148  (*B3) [i] = (A3->size() + i) * 100;
149  }
150 
151  for ( UInt i (0) ; i < static_cast<UInt> (A4->size() ) ; ++i )
152  {
153  (*A4) [i] = i * 10;
154  }
155  for ( UInt i (0) ; i < static_cast<UInt> (B4->size() ) ; ++i )
156  {
157  (*B4) [i] = (A4->size() + i) * 1000;
158  }
159 
160  UInt result = TestFunction ( A1, B1, A2, B2, A3, B3, A4, B4 );
161 
162 
163 
164  // ===================================================
165  // TEST: BOOST::NUMERIC::UBLAS::VECTOR
166  // ===================================================
167 
168  typedef boost::numeric::ublas::vector<Real> VectorBOOST;
169  typedef std::shared_ptr<VectorBOOST> VectorBOOST_ptr;
170 
171  // BASE Vectors
172  VectorBOOST_ptr A1_BOOST, B1_BOOST, A2_BOOST, B2_BOOST;
173 
174  A1_BOOST.reset ( new VectorBOOST (2) );
175  (*A1_BOOST) [0] = 0;
176  (*A1_BOOST) [1] = 1;
177 
178  B1_BOOST.reset ( new VectorBOOST (3) );
179  (*B1_BOOST) [0] = 2;
180  (*B1_BOOST) [1] = 3;
181  (*B1_BOOST) [2] = 4;
182 
183  A2_BOOST.reset ( new VectorBOOST (2) );
184  (*A2_BOOST) [0] = 0;
185  (*A2_BOOST) [1] = 10;
186 
187  B2_BOOST.reset ( new VectorBOOST (3) );
188  (*B2_BOOST) [0] = 20;
189  (*B2_BOOST) [1] = 30;
190  (*B2_BOOST) [2] = 40;
191 
192  // NOTE:
193  // To perform a test with BOOST we first have to implement a derived class
194  // from VectorContainer in which re-implement operator= and operator*
195 
196  //TestFunction( A1_BOOST, B1_BOOST, A2_BOOST, B2_BOOST );
197 
198 #ifdef HAVE_MPI
199  std::cout << std::endl << "MPI Finalization" << std::endl;
200  MPI_Finalize();
201 #endif
202 
203  return ( result );
204 }
VectorEpetra - The Epetra Vector format Wrapper.
int32_type Int
Generic integer data.
Definition: LifeV.hpp:188
void updateInverseJacobian(const UInt &iQuadPt)
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
int main(int argc, char **argv)
Definition: dummy.cpp:5
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191