LifeV
ETA_ADR2DTest.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_ADR2DTest.cpp
28  @author L. Pasquale <luca.pasquale@mail.polimi.it>
29  @date 2012-11-20
30  */
31 
32 // ===================================================
33 //! Includes
34 // ===================================================
35 
36 #include "ETA_ADR2DTest.hpp"
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 Real fRhs ( const Real& /* t */, const Real& /*x*/, const Real& /*y*/, const Real& /* z */ , const ID& /* i */ )
78 {
79  return 2;
80 }
81 
82 
83 // ===================================================
84 //! Constructors
85 // ===================================================
86 
87 ETA_ADR2DTest::ETA_ADR2DTest ()
88 {
89 
90 #ifdef EPETRA_MPI
91  M_comm.reset ( new Epetra_MpiComm ( MPI_COMM_WORLD ) );
92 #else
93  M_comm.reset ( new Epetra_SerialComm() );
94 #endif
95 
96 }
97 
98 // ===================================================
99 //! Methods
100 // ===================================================
101 
102 Real
103 ETA_ADR2DTest::run()
104 {
105  bool verbose (M_comm->MyPID() == 0);
106  // ---------------------------------------------------------------
107  // We define the mesh and parition it. We use the domain
108  // (-1,1)x(-1,1) and a structured mesh.
109  // ---------------------------------------------------------------
110 
111  if (verbose)
112  {
113  std::cout << " -- Building and partitioning the mesh ... " << std::flush;
114  }
115 
116  const UInt Nelements (10);
117 
118  std::shared_ptr< mesh_Type > fullMeshPtr (new mesh_Type);
119 
120  regularMesh2D ( *fullMeshPtr, 0, Nelements, Nelements, false,
121  2.0, 2.0,
122  -0.0, -0.0);
123 
124  std::shared_ptr< mesh_Type > meshPtr;
125  {
126  MeshPartitioner< mesh_Type > meshPart (fullMeshPtr, M_comm);
127  meshPtr = meshPart.meshPartition();
128  }
129 
130  fullMeshPtr.reset();
131 
132  if (verbose)
133  {
134  std::cout << " done ! " << std::endl;
135  }
136 
137 
138  // ---------------------------------------------------------------
139  // We start by defining the finite element spaces. We use the type
140  // FESpace for the classical way and copy it in an ETFESpace for
141  // the ET assembly.
142  //
143  // In the end, both solution spaces display their respective
144  // number of degrees of freedom, which must be the same.
145  // ---------------------------------------------------------------
146 
147  if (verbose)
148  {
149  std::cout << " -- Building FESpaces ... " << std::flush;
150  }
151 
152  std::string uOrder ("P1");
153  std::string bOrder ("P1");
154 
155  std::shared_ptr<FESpace< mesh_Type, MapEpetra > > uSpace
156  ( new FESpace< mesh_Type, MapEpetra > (meshPtr, uOrder, 1, M_comm) );
157 
158  std::shared_ptr<FESpace< mesh_Type, MapEpetra > > betaSpace
159  ( new FESpace< mesh_Type, MapEpetra > (meshPtr, bOrder, 2, M_comm) );
160 
161  if (verbose)
162  {
163  std::cout << " done ! " << std::endl;
164  }
165  if (verbose)
166  {
167  std::cout << " ---> Dofs: " << uSpace->dof().numTotalDof() << std::endl;
168  }
169 
170  if (verbose)
171  {
172  std::cout << " -- Building ETFESpaces ... " << std::flush;
173  }
174 
175  std::shared_ptr<ETFESpace< mesh_Type, MapEpetra, 2, 1 > > ETuSpace
176  ( new ETFESpace< mesh_Type, MapEpetra, 2, 1 > (meshPtr, & (uSpace->refFE() ), & (uSpace->fe().geoMap() ), M_comm) );
177 
178  std::shared_ptr<ETFESpace< mesh_Type, MapEpetra, 2, 2 > > ETbetaSpace
179  ( new ETFESpace< mesh_Type, MapEpetra, 2, 2 > (meshPtr, & (betaSpace->refFE() ), & (betaSpace->fe().geoMap() ), M_comm) );
180 
181  if (verbose)
182  {
183  std::cout << " done ! " << std::endl;
184  }
185  if (verbose)
186  {
187  std::cout << " ---> Dofs: " << ETuSpace->dof().numTotalDof() << std::endl;
188  }
189 
190 
191  // ---------------------------------------------------------------
192  // We then interpolate the advection function.
193  // This is can only be performed with the classical FESpace.
194  // ---------------------------------------------------------------
195 
196  if (verbose)
197  {
198  std::cout << " -- Interpolating the advection field ... " << std::flush;
199  }
200 
201  vector_Type beta (betaSpace->map(), Repeated);
202  betaSpace->interpolate (static_cast<FESpace< mesh_Type, MapEpetra >::function_Type> (betaFct), beta, 0.0);
203 
204  if (verbose)
205  {
206  std::cout << " done! " << std::endl;
207  }
208 
209 
210  // ---------------------------------------------------------------
211  // We build now two matrices, one for each assembly procedure in
212  // order to compare them in the end. There is no difference in the
213  // construction of the matrices.
214  // ---------------------------------------------------------------
215 
216  if (verbose)
217  {
218  std::cout << " -- Defining the matrices ... " << std::flush;
219  }
220 
221  std::shared_ptr<matrix_Type> systemMatrix (new matrix_Type ( uSpace->map() ) );
222  *systemMatrix *= 0.0;
223 
224  std::shared_ptr<matrix_Type> ETsystemMatrix (new matrix_Type ( ETuSpace->map() ) );
225  *ETsystemMatrix *= 0.0;
226 
227  if (verbose)
228  {
229  std::cout << " done! " << std::endl;
230  }
231 
232  // ---------------------------------------------------------------
233  // Definition of the RHS
234  // ---------------------------------------------------------------
235 
236  if (verbose)
237  {
238  std::cout << " -- Defining and interpolating the RHS ... " << std::flush;
239  }
240 
241  vector_Type rhs (uSpace->map(), Repeated);
242  rhs *= 0.0;
243  vector_Type ETrhs (uSpace->map(), Repeated);
244  ETrhs *= 0.0;
245 
246  vector_Type fInterpolated (uSpace->map(), Repeated);
247  fInterpolated *= 0.0;
248  uSpace->interpolate ( static_cast<FESpace< mesh_Type, MapEpetra >::function_Type> (fRhs) , fInterpolated, 0.0 );
249 
250  if (verbose)
251  {
252  std::cout << " done! " << std::endl;
253  }
254 
255  // ---------------------------------------------------------------
256  // We come now to the assembly itself. We start with the classical
257  // way, which consists in instentiating an ADRAssembler and call
258  // the different methods from that class.
259  //
260  // Using a LifeChrono, we monitor the time required for the
261  // assembly.
262  //
263  // The terms assembled correspond, in the order of appearance, to
264  // - a laplacian term
265  // - an advection term
266  // - a reaction term (with coefficient 2.0), aka mass term.
267  //
268  // We also assemble the RHS
269  // ---------------------------------------------------------------
270 
271  LifeChrono StdChrono;
272  StdChrono.start();
273 
274  if (verbose)
275  {
276  std::cout << " -- Classical assembly ... " << std::flush;
277  }
278 
280 
281  adrAssembler.setup (uSpace, betaSpace);
282 
283  adrAssembler.addDiffusion (systemMatrix, 1.0);
284 
285  adrAssembler.addAdvection (systemMatrix, beta);
286 
287  adrAssembler.addMass (systemMatrix, 2.0);
288 
289  adrAssembler.addMassRhs (rhs, fInterpolated);
290 
291  StdChrono.stop();
292 
293  if (verbose)
294  {
295  std::cout << " done ! " << std::endl;
296  }
297  if (verbose)
298  {
299  std::cout << " Time : " << StdChrono.diff() << std::endl;
300  }
301 
302 
303  // ---------------------------------------------------------------
304  // We perform now the same assembly with the ET assembly and still
305  // monitor the timings required.
306  //
307  // Remark also that the Real constants and VectorSmall constants
308  // can be used directly in the expressions to integrate.
309  // ---------------------------------------------------------------
310 
311  LifeChrono ETChrono;
312  ETChrono.start();
313 
314  if (verbose)
315  {
316  std::cout << " -- ET assembly ... " << std::flush;
317  }
318 
319 
320  {
321  using namespace ExpressionAssembly;
322 
323  integrate ( elements (ETuSpace->mesh() ),
324  uSpace->qr(),
325  ETuSpace,
326  ETuSpace,
327 
328  dot ( grad (phi_i) , grad (phi_j) )
329  + dot ( grad (phi_j) , value (ETbetaSpace, beta) ) *phi_i
330  + 2.0 * phi_i * phi_j
331 
332  )
333  >> ETsystemMatrix;
334 
335  integrate ( elements (ETuSpace->mesh() ),
336  uSpace->qr(),
337  ETuSpace,
338 
339  value (ETuSpace, fInterpolated) *phi_i
340 
341  )
342  >> ETrhs;
343 
344  }
345 
346  ETChrono.stop();
347 
348  if (verbose)
349  {
350  std::cout << " done! " << std::endl;
351  }
352  if (verbose)
353  {
354  std::cout << " Time : " << ETChrono.diff() << std::endl;
355  }
356 
357 
358  // ---------------------------------------------------------------
359  // We finally need to check that both yield the same matrix. In
360  // that aim, we need to finalize both matrices.
361  // ---------------------------------------------------------------
362 
363  if (verbose)
364  {
365  std::cout << " -- Closing the matrices and vectors... " << std::flush;
366  }
367 
368  systemMatrix->globalAssemble();
369  ETsystemMatrix->globalAssemble();
370 
371  rhs.globalAssemble();
372  ETrhs.globalAssemble();
373 
374  if (verbose)
375  {
376  std::cout << " done ! " << std::endl;
377  }
378 
379  // ---------------------------------------------------------------
380  // We compute now the matrix of the difference and finally the
381  // norm of the difference. This should be very low if the two
382  // matrices are identical.
383  // ---------------------------------------------------------------
384 
385  if (verbose)
386  {
387  std::cout << " -- Computing the error ... " << std::flush;
388  }
389 
390  std::shared_ptr<matrix_Type> checkMatrix (new matrix_Type ( ETuSpace->map() ) );
391  *checkMatrix *= 0.0;
392 
393  *checkMatrix += *systemMatrix;
394  *checkMatrix += (*ETsystemMatrix) * (-1);
395 
396  checkMatrix->globalAssemble();
397 
398  Real errorNorm ( checkMatrix->normInf() );
399 
400 
401  // ---------------------------------------------------------------
402  // Finally, we check that the two right hand sides are similar
403  // with the two assembly procedures. We also check that the two
404  // integrals correspond to their exact values.
405  // ---------------------------------------------------------------
406 
407  vector_Type checkRhs (ETuSpace->map(), Repeated);
408  checkRhs = 0.0;
409 
410  checkRhs += rhs;
411  checkRhs -= ETrhs;
412 
413  checkRhs.globalAssemble();
414 
415  vector_Type checkRhsUnique (checkRhs, Unique);
416 
417  Real errorRhsNorm ( checkRhsUnique.normInf() );
418 
419  if (verbose)
420  {
421  std::cout << " done ! " << std::endl;
422  }
423 
424  if (verbose)
425  {
426  std::cout << " Matrix Error : " << errorNorm << std::endl;
427  }
428  if (verbose)
429  {
430  std::cout << " Rhs error : " << errorRhsNorm << std::endl;
431  }
432 
433 
434  // ---------------------------------------------------------------
435  // Finally we return the maximum of the two errors
436  // ---------------------------------------------------------------
437 
438 
439  return std::max (errorNorm, errorRhsNorm);
440 
441 } // run
Real fRhs(const Real &, const Real &, const Real &, const Real &, const ID &i)
RegionMesh< LinearTriangle > mesh_Type
uint32_type ID
IDs.
Definition: LifeV.hpp:194
MatrixEpetra< Real > matrix_Type
VectorEpetra vector_Type
double Real
Generic real data.
Definition: LifeV.hpp:175
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.