LifeV
lifev/eta/testsuite/mt_assembly/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 Test for building matrices with a static graph and ETA
30 
31  @author Radu Popescu <radu.popescu@epfl.ch>
32  @date 2012-03-19
33  */
34 
35 #pragma GCC diagnostic ignored "-Wunused-variable"
36 #pragma GCC diagnostic ignored "-Wunused-parameter"
37 
38 #include <cstdlib>
39 
40 #include <Epetra_ConfigDefs.h>
41 #ifdef EPETRA_MPI
42 #include <mpi.h>
43 #include <Epetra_MpiComm.h>
44 #else
45 #include <Epetra_SerialComm.h>
46 #endif
47 
48 #include <Epetra_FECrsGraph.h>
49 
50 #pragma GCC diagnostic warning "-Wunused-variable"
51 #pragma GCC diagnostic warning "-Wunused-parameter"
52 
53 #include <lifev/core/LifeV.hpp>
54 #include <lifev/core/util/WallClock.hpp>
55 #include <lifev/core/util/OpenMPParameters.hpp>
56 
57 #include <lifev/core/mesh/MeshPartitioner.hpp>
58 #include <lifev/core/mesh/RegionMesh3DStructured.hpp>
59 #include <lifev/core/mesh/RegionMesh.hpp>
60 
61 #include <lifev/core/array/MatrixEpetra.hpp>
62 
63 #include <lifev/eta/fem/ETFESpace.hpp>
64 
65 #include <lifev/eta/expression/Integrate.hpp>
66 #include <lifev/eta/expression/BuildGraph.hpp>
67 
68 
69 #include <boost/shared_ptr.hpp>
70 
71 
72 using namespace LifeV;
73 
75 typedef MatrixEpetra<Real> matrix_Type;
76 
77 int main ( int argc, char** argv )
78 {
79 
80 #ifdef HAVE_MPI
81  MPI_Init (&argc, &argv);
82  std::shared_ptr<Epetra_Comm> Comm (new Epetra_MpiComm (MPI_COMM_WORLD) );
83 #else
84  std::shared_ptr<Epetra_Comm> Comm (new Epetra_SerialComm);
85 #endif
86 
87  const bool verbose (Comm->MyPID() == 0);
88 
89  if (argc != 5)
90  {
91  if (verbose)
92  {
93  std::cout << "Please run program as " << argv[0]
94  << " " << "<num_elements> " << "<num_threads>"
95  << " " << "<scheduler> " << "<chunkSize>\n";
96  return EXIT_FAILURE;
97  }
98  }
99 
100  OpenMPParameters ompParams;
101 
102  const UInt Nelements = std::atoi (argv[1]);
103  ompParams.numThreads = std::atoi (argv[2]);
104 #ifdef _OPENMP
105  ompParams.scheduler = static_cast<omp_sched_t> (std::atoi (argv[3]) );
106 #endif
107  ompParams.chunkSize = std::atoi (argv[4]);
108 
109  if (verbose)
110  {
111  std::cout << " -- Building and partitioning the mesh ... " << std::flush;
112  }
113 
114 
115  std::shared_ptr< mesh_Type > fullMeshPtr (new mesh_Type);
116 
117  regularMesh3D ( *fullMeshPtr, 1, Nelements, Nelements, Nelements, false,
118  2.0, 2.0, 2.0,
119  -1.0, -1.0, -1.0);
120 
121  // MeshPartitioner< mesh_Type > meshPart (fullMeshPtr, Comm);
122  MeshPartitioner< mesh_Type > meshPart;
123  meshPart.setPartitionOverlap ( 1 );
124  meshPart.doPartition ( fullMeshPtr, Comm );
125 
126  fullMeshPtr.reset();
127 
128  if (verbose)
129  {
130  std::cout << " done ! " << std::endl;
131  }
132 
133 
134  if (verbose)
135  {
136  std::cout << " -- Building ETFESpaces ... " << std::flush;
137  }
138 
139  std::shared_ptr<ETFESpace< mesh_Type, MapEpetra, 3, 1 > > uSpace
140  ( new ETFESpace< mesh_Type, MapEpetra, 3, 1 > (meshPart, &feTetraP1, Comm) );
141 
142  if (verbose)
143  {
144  std::cout << " done ! " << std::endl;
145  }
146  if (verbose)
147  {
148  std::cout << " ---> Dofs: " << uSpace->dof().numTotalDof() << std::endl;
149  }
150 
151  std::shared_ptr<matrix_Type> closedSystemMatrix;
152 
153  if (verbose)
154  {
155  std::cout << " -- Precomputing matrix graph ... " << std::flush;
156  }
157 
158  WallClock timer;
159  timer.start();
160  std::shared_ptr<Epetra_FECrsGraph> matrixGraph;
161  {
162  using namespace ExpressionAssembly;
163 
164  // We first build a static graph for our problem matrix
165  matrixGraph.reset (new Epetra_FECrsGraph (Copy, * (uSpace->map().map (Unique) ), 0, true) );
166 
167  buildGraph ( elements (uSpace->mesh() ),
168  quadRuleTetra4pt,
169  uSpace,
170  uSpace,
171  dot ( grad (phi_i) , grad (phi_j) )
172  ) >> matrixGraph;
173  }
174  matrixGraph->GlobalAssemble();
175  timer.stop();
176 
177  if (verbose)
178  {
179  std::cout << " done in " << timer.elapsedTime() << "s." << std::endl;
180  }
181 
182  if (verbose)
183  {
184  std::cout << " -- Assembling the Laplace matrix with a precomputed graph ... " << std::flush;
185  }
186 
187  timer.start();
188  {
189  using namespace ExpressionAssembly;
190 
191  // We build the system matrix using the precomputed graph
192  closedSystemMatrix.reset (new matrix_Type ( uSpace->map(), *matrixGraph , true) );
193  *closedSystemMatrix *= 0.0;
194 
195  // Finally, we perform the FE assembly, which should be faster with
196  // a closed matrix
197  integrate ( elements (uSpace->mesh() ),
198  quadRuleTetra4pt,
199  uSpace,
200  uSpace,
201  dot ( grad (phi_i) , grad (phi_j) ), ompParams
202  ) >> closedSystemMatrix;
203 
204  closedSystemMatrix->globalAssemble();
205  }
206  timer.stop();
207 
208  if (verbose)
209  {
210  std::cout << " done in " << timer.elapsedTime() << "s." << std::endl;
211  }
212 
213  Real closedMatrixNorm ( closedSystemMatrix->normInf() );
214 
215  if (verbose)
216  {
217  std::cout << " Closed matrix norm : " << closedMatrixNorm << std::endl;
218  }
219 
220 #ifdef HAVE_MPI
221  MPI_Finalize();
222 #endif
223 
224  Real closedMatrixNormDiff (std::abs (closedMatrixNorm - 3.2) );
225 
226  if (verbose)
227  {
228  std::cout << " Error (closed): " << closedMatrixNormDiff << std::endl;
229  }
230 
231  Real testTolerance (1e-10);
232 
233  if ( closedMatrixNormDiff >= testTolerance )
234  {
235  return ( EXIT_FAILURE );
236  }
237  return ( EXIT_SUCCESS );
238 
239 }
int main(int argc, char **argv)
MatrixEpetra< Real > matrix_Type
double Real
Generic real data.
Definition: LifeV.hpp:175
RegionMesh< LinearTetra > mesh_Type