LifeV
ETA_VectorialADR2DTest.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  @file ETA_VectorialADR2DTest.cpp
28  @author L. Pasquale <luca.pasquale@mail.polimi.it>
29  @date 2012-11-20
30  */
31 
32 // ===================================================
33 //! Includes
34 // ===================================================
35 
37 
38 #include <lifev/eta/fem/ETFESpace.hpp>
39 #include <lifev/eta/expression/Integrate.hpp>
40 
41 
42 // ===================================================
43 //! Namespaces & define
44 // ===================================================
45 
46 // ---------------------------------------------------------------
47 // We work in the LifeV namespace and define the mesh, matrix and
48 // vector types that we will need several times.
49 // ---------------------------------------------------------------
50 
51 using namespace LifeV;
52 
54 typedef MatrixEpetra<Real> matrix_Type;
55 typedef VectorEpetra vector_Type;
56 
57 // ---------------------------------------------------------------
58 // We define then a function, which represents a velocity field,
59 // of magnitude 1 in the y direction (supposing an (x,y)
60 // reference frame).
61 // ---------------------------------------------------------------
62 
63 
64 // ===================================================
65 //! Functions
66 // ===================================================
67 
68 Real betaFct ( const Real& /* t */, const Real& /* x */, const Real& /* y */, const Real& /* z */, const ID& i )
69 {
70  if (i == 1)
71  {
72  return 1;
73  }
74  return 0;
75 }
76 
77 // ===================================================
78 //! Constructors
79 // ===================================================
80 
81 ETA_VectorialADR2DTest::ETA_VectorialADR2DTest ()
82 {
83 
84 #ifdef EPETRA_MPI
85  M_comm.reset ( new Epetra_MpiComm ( MPI_COMM_WORLD ) );
86 #else
87  M_comm.reset ( new Epetra_SerialComm() );
88 #endif
89 
90 }
91 
92 // ===================================================
93 //! Methods
94 // ===================================================
95 
96 Real
97 ETA_VectorialADR2DTest::run()
98 {
99  bool verbose (M_comm->MyPID() == 0);
100  // ---------------------------------------------------------------
101  // We define the mesh and parition it. We use the domain
102  // (-1,1)x(-1,1) and a structured mesh.
103  // ---------------------------------------------------------------
104 
105  if (verbose)
106  {
107  std::cout << " -- Building and partitioning the mesh ... " << std::flush;
108  }
109 
110  const UInt Nelements (10);
111 
112  std::shared_ptr< mesh_Type > fullMeshPtr (new mesh_Type);
113 
114  regularMesh2D ( *fullMeshPtr, 0, Nelements, Nelements, false,
115  2.0, 2.0,
116  -0.0, -0.0);
117 
118  std::shared_ptr< mesh_Type > meshPtr;
119  {
120  MeshPartitioner< mesh_Type > meshPart (fullMeshPtr, M_comm);
121  meshPtr = meshPart.meshPartition();
122  }
123 
124  fullMeshPtr.reset();
125 
126  if (verbose)
127  {
128  std::cout << " done ! " << std::endl;
129  }
130 
131 
132  // ---------------------------------------------------------------
133  // We start by defining the finite element spaces. We use the type
134  // FESpace for the classical way and copy it in an ETFESpace for
135  // the ET assembly.
136  //
137  // In the end, both solution spaces display their respective
138  // number of degrees of freedom, which must be the same.
139  // ---------------------------------------------------------------
140 
141  if (verbose)
142  {
143  std::cout << " -- Building FESpaces ... " << std::flush;
144  }
145 
146  std::string uOrder ("P1");
147  std::string bOrder ("P1");
148 
149  std::shared_ptr<FESpace< mesh_Type, MapEpetra > > uSpace
150  ( new FESpace< mesh_Type, MapEpetra > (meshPtr, uOrder, 2, M_comm) );
151 
152  std::shared_ptr<FESpace< mesh_Type, MapEpetra > > betaSpace
153  ( new FESpace< mesh_Type, MapEpetra > (meshPtr, bOrder, 2, M_comm) );
154 
155  if (verbose)
156  {
157  std::cout << " done ! " << std::endl;
158  }
159  if (verbose)
160  {
161  std::cout << " ---> Dofs: " << uSpace->dof().numTotalDof() << std::endl;
162  }
163 
164  if (verbose)
165  {
166  std::cout << " -- Building ETFESpaces ... " << std::flush;
167  }
168 
169  std::shared_ptr<ETFESpace< mesh_Type, MapEpetra, 2, 2 > > ETuSpace
170  ( new ETFESpace< mesh_Type, MapEpetra, 2, 2 > (meshPtr, & (uSpace->refFE() ), & (uSpace->fe().geoMap() ), M_comm) );
171 
172  std::shared_ptr<ETFESpace< mesh_Type, MapEpetra, 2, 2 > > ETbetaSpace
173  ( new ETFESpace< mesh_Type, MapEpetra, 2, 2 > (meshPtr, & (betaSpace->refFE() ), & (betaSpace->fe().geoMap() ), M_comm) );
174 
175  if (verbose)
176  {
177  std::cout << " done ! " << std::endl;
178  }
179  if (verbose)
180  {
181  std::cout << " ---> Dofs: " << ETuSpace->dof().numTotalDof() << std::endl;
182  }
183 
184 
185  // ---------------------------------------------------------------
186  // We interpolate then the advection function of the mesh at hand.
187  // This is performed with the classical FESpace only.
188  // ---------------------------------------------------------------
189 
190  if (verbose)
191  {
192  std::cout << " -- Interpolating the advection field ... " << std::flush;
193  }
194 
195  vector_Type beta (betaSpace->map(), Repeated);
196  betaSpace->interpolate (static_cast<FESpace< mesh_Type, MapEpetra >::function_Type> (betaFct), beta, 0.0);
197 
198  if (verbose)
199  {
200  std::cout << " done! " << std::endl;
201  }
202 
203 
204  // ---------------------------------------------------------------
205  // We build now two matrices, one for each assembly procedure in
206  // order to compare them in the end. There is no difference in the
207  // construction of the matrices.
208  // ---------------------------------------------------------------
209 
210  if (verbose)
211  {
212  std::cout << " -- Defining the matrices ... " << std::flush;
213  }
214 
215  std::shared_ptr<matrix_Type> systemMatrix (new matrix_Type ( uSpace->map() ) );
216  *systemMatrix *= 0.0;
217 
218  std::shared_ptr<matrix_Type> ETsystemMatrix (new matrix_Type ( ETuSpace->map() ) );
219  *ETsystemMatrix *= 0.0;
220 
221  if (verbose)
222  {
223  std::cout << " done! " << std::endl;
224  }
225 
226 
227  // ---------------------------------------------------------------
228  // We come now to the assembly itself. We start with the classical
229  // way, which consists in instentiating an ADRAssembler and call
230  // the different methods from that class.
231  //
232  // Using a LifeChrono, we monitor the time required for the
233  // assembly.
234  //
235  // The terms assembled correspond, in the order of appearance, to
236  // - a laplacian term
237  // - an advection term
238  // - a reaction term (with coefficient 2.0), aka mass term.
239  // ---------------------------------------------------------------
240 
241  LifeChrono StdChrono;
242  StdChrono.start();
243 
244  if (verbose)
245  {
246  std::cout << " -- Classical assembly ... " << std::flush;
247  }
248 
250 
251  adrAssembler.setup (uSpace, betaSpace);
252 
253  adrAssembler.addDiffusion (systemMatrix, 1.0);
254 
255  adrAssembler.addAdvection (systemMatrix, beta);
256 
257  adrAssembler.addMass (systemMatrix, 2.0);
258 
259  StdChrono.stop();
260 
261  if (verbose)
262  {
263  std::cout << " done ! " << std::endl;
264  }
265  if (verbose)
266  {
267  std::cout << " Time : " << StdChrono.diff() << std::endl;
268  }
269 
270 
271  // ---------------------------------------------------------------
272  // We perform now the same assembly with the ET assembly and still
273  // monitor the timings required.
274  //
275  // Remark also that the Real constants and VectorSmall constants
276  // can be used directly in the expressions to integrate.
277  // ---------------------------------------------------------------
278 
279  LifeChrono ETChrono;
280  ETChrono.start();
281 
282  if (verbose)
283  {
284  std::cout << " -- ET assembly ... " << std::flush;
285  }
286 
287 
288  {
289  using namespace ExpressionAssembly;
290 
291  integrate ( elements (ETuSpace->mesh() ),
292  uSpace->qr(),
293  ETuSpace,
294  ETuSpace,
295 
296  dot ( grad (phi_i) , grad (phi_j) )
297  + dot ( grad (phi_j) * value (ETbetaSpace, beta), phi_i)
298  + 2.0 * dot (phi_i, phi_j)
299 
300  )
301  >> ETsystemMatrix;
302 
303  }
304 
305  ETChrono.stop();
306 
307  if (verbose)
308  {
309  std::cout << " done! " << std::endl;
310  }
311  if (verbose)
312  {
313  std::cout << " Time : " << ETChrono.diff() << std::endl;
314  }
315 
316 
317  // ---------------------------------------------------------------
318  // We finally need to check that both yield the same matrix. In
319  // that aim, we need to finalize both matrices.
320  // ---------------------------------------------------------------
321 
322  if (verbose)
323  {
324  std::cout << " -- Closing the matrices ... " << std::flush;
325  }
326 
327  systemMatrix->globalAssemble();
328  ETsystemMatrix->globalAssemble();
329 
330  if (verbose)
331  {
332  std::cout << " done ! " << std::endl;
333  }
334 
335  // ---------------------------------------------------------------
336  // We compute now the matrix of the difference and finally the
337  // norm of the difference. This should be very low if the two
338  // matrices are identical.
339  // ---------------------------------------------------------------
340 
341  if (verbose)
342  {
343  std::cout << " -- Computing the error ... " << std::flush;
344  }
345 
346  std::shared_ptr<matrix_Type> checkMatrix (new matrix_Type ( ETuSpace->map() ) );
347  *checkMatrix *= 0.0;
348 
349  *checkMatrix += *systemMatrix;
350  *checkMatrix += (*ETsystemMatrix) * (-1);
351 
352  checkMatrix->globalAssemble();
353 
354  Real errorNorm ( checkMatrix->normInf() );
355 
356  if (verbose)
357  {
358  std::cout << " done ! " << std::endl;
359  }
360 
361 
362 
363  // ---------------------------------------------------------------
364  // We finally display the error norm and compare with the
365  // tolerance of the test.
366  // ---------------------------------------------------------------
367 
368  if (verbose)
369  {
370  std::cout << " Matrix Error : " << errorNorm << std::endl;
371  }
372 
373  return errorNorm;
374 
375 } // run
uint32_type ID
IDs.
Definition: LifeV.hpp:194
MatrixEpetra< Real > matrix_Type
double Real
Generic real data.
Definition: LifeV.hpp:175
VectorEpetra vector_Type
RegionMesh< LinearTriangle > mesh_Type
ADRAssembler - This class add into given matrices terms corresponding to the space discretization of ...
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191
Real betaFct(const Real &, const Real &, const Real &, const Real &, const ID &i)
Functions.