LifeV
lifev/eta/tutorials/12_ETA_LaplacianPhiJ/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 Tutorial for laplacian of trial functions.
30 
31  @author Davide Forti <davide.forti@epfl.ch>
32  @date 2013-05
33 
34  Test to check (and to show) how the laplacian pf the test
35  and trial functions is computed
36 
37  */
38 
39 // ---------------------------------------------------------------
40 // We reuse the same files as in the previous tutorial. We add the
41 // chrono to measure the timings.
42 // ---------------------------------------------------------------
43 
44 
45 #include <Epetra_ConfigDefs.h>
46 #ifdef EPETRA_MPI
47 #include <mpi.h>
48 #include <Epetra_MpiComm.h>
49 #else
50 #include <Epetra_SerialComm.h>
51 #endif
52 
53 
54 #include <lifev/core/LifeV.hpp>
55 
56 #include <lifev/core/mesh/MeshPartitioner.hpp>
57 #include <lifev/core/mesh/RegionMesh3DStructured.hpp>
58 #include <lifev/core/mesh/RegionMesh.hpp>
59 #include <lifev/core/mesh/MeshData.hpp>
60 
61 #include <lifev/core/array/MatrixEpetra.hpp>
62 #include <lifev/core/array/VectorEpetra.hpp>
63 
64 #include <lifev/eta/fem/ETFESpace.hpp>
65 #include <lifev/eta/expression/Integrate.hpp>
66 
67 #include <lifev/core/fem/FESpace.hpp>
68 
69 #include <lifev/core/util/LifeChrono.hpp>
70 
71 using namespace LifeV;
72 
74 typedef MatrixEpetra<Real> matrix_Type;
75 typedef VectorEpetra vector_Type;
76 
77 Real uTestFunctions ( const Real& /* t */, const Real& x , const Real& y , const Real& z , const ID& i )
78 {
79  if ( i == 0 )
80  return x*x*0.5;
81  else if ( i == 1)
82  return y*y*0.5;
83  else
84  return z*z*0.5;
85 }
86 
87 int main ( int argc, char** argv )
88 {
89 
90 #ifdef HAVE_MPI
91  MPI_Init (&argc, &argv);
92  std::shared_ptr<Epetra_Comm> Comm (new Epetra_MpiComm (MPI_COMM_WORLD) );
93 #else
94  std::shared_ptr<Epetra_Comm> Comm (new Epetra_SerialComm);
95 #endif
96 
97  const bool verbose (Comm->MyPID() == 0);
98 
99  const UInt Nelements (10);
100 
101  std::shared_ptr< mesh_Type > fullMeshPtr (new mesh_Type);
102 
103  Real length = 3.0;
104 
105 
106  regularMesh3D ( *fullMeshPtr, 1, Nelements, Nelements, Nelements, false,
107  length, length, length,
108  0.0, 0.0, 0.0);
109 
110  MeshPartitioner< mesh_Type > meshPart (fullMeshPtr, Comm);
111  std::shared_ptr< mesh_Type > meshPtr (meshPart.meshPartition() );
112 
113  fullMeshPtr.reset();
114 
115  std::string uOrder ("P2");
116 
117  ///////////////////////////////////
118  // Testing the scalar field case //
119  ///////////////////////////////////
120 
121  std::shared_ptr<FESpace< mesh_Type, MapEpetra > > uSpace
122  ( new FESpace< mesh_Type, MapEpetra > (meshPtr, uOrder, 1, Comm) );
123 
124  std::shared_ptr<ETFESpace< mesh_Type, MapEpetra, 3, 1 > > ETuSpace
125  ( new ETFESpace< mesh_Type, MapEpetra, 3, 1 > (meshPtr, & (uSpace->refFE() ), & (uSpace->fe().geoMap() ), Comm) );
126 
127  vector_Type vectorTestFunctions (uSpace->map(), Unique);
128  uSpace->interpolate (static_cast<FESpace< mesh_Type, MapEpetra >::function_Type> (uTestFunctions), vectorTestFunctions, 0.0);
129 
130  matrix_Type integral (uSpace->map() );
131 
132  {
133  using namespace ExpressionAssembly;
134 
135  integrate ( elements (ETuSpace->mesh() ),
136  uSpace->qr(),
137  ETuSpace,
138  ETuSpace,
139  laplacian(phi_i) * laplacian(phi_j)
140  )
141  >> integral;
142  }
143 
144  integral.globalAssemble();
145 
146  Real result = 0.0;
147 
148  vector_Type first(uSpace->map());
149  first = integral * vectorTestFunctions;
150  result = vectorTestFunctions.dot(first);
151 
152  if ( Comm->MyPID() == 0 )
153  {
154  std::cout << "\n\nSCALAR CASE " << std::endl;
155  std::cout << "\nThe volume is = " << 27 << std::endl;
156  std::cout << "\nThe result is = " << result << std::endl;
157  std::cout << "\nThe error is = " << result-27 << std::endl;
158  }
159 
160  Real error_scalar = result-27.0;
161 
162  ///////////////////////////////////
163  // Testing the vector field case //
164  ///////////////////////////////////
165 
166  std::shared_ptr<FESpace< mesh_Type, MapEpetra > > uSpaceVec
167  ( new FESpace< mesh_Type, MapEpetra > (meshPtr, uOrder, 3, Comm) );
168 
169  std::shared_ptr<ETFESpace< mesh_Type, MapEpetra, 3, 3 > > ETuSpaceVec
170  ( new ETFESpace< mesh_Type, MapEpetra, 3, 3 > (meshPtr, & (uSpaceVec->refFE() ), & (uSpaceVec->fe().geoMap() ), Comm) );
171 
172  vector_Type vectorTestFunctionsVec (uSpaceVec->map(), Unique);
173  uSpaceVec->interpolate (static_cast<FESpace< mesh_Type, MapEpetra >::function_Type> (uTestFunctions), vectorTestFunctionsVec, 0.0);
174 
175  matrix_Type integralVec (uSpaceVec->map() );
176 
177  {
178  using namespace ExpressionAssembly;
179 
180  integrate ( elements (ETuSpaceVec->mesh() ),
181  uSpaceVec->qr(),
182  ETuSpaceVec,
183  ETuSpaceVec,
184  dot ( laplacian(phi_i), laplacian(phi_j) )
185  )
186  >> integralVec;
187  }
188 
189  integralVec.globalAssemble();
190 
191  result = 0.0;
192 
193  vector_Type firstVec(uSpaceVec->map());
194  firstVec = integralVec * vectorTestFunctionsVec;
195  result = vectorTestFunctionsVec.dot(firstVec);
196 
197  if ( Comm->MyPID() == 0 )
198  {
199  std::cout << "\n\nVECTORIAL CASE " << std::endl;
200  std::cout << "\nThe volume is = " << 81 << std::endl;
201  std::cout << "\nThe result is = " << result << std::endl;
202  std::cout << "\nThe error is = " << result-81 << "\n\n";
203  }
204 
205  Real error_vectorial = result-81.0;
206 
207 #ifdef HAVE_MPI
208  MPI_Finalize();
209 #endif
210 
211  if ( std::fabs(error_scalar) < 1e-9 && std::fabs(error_vectorial) < 1e-9 )
212  {
213  return ( EXIT_SUCCESS );
214  }
215 
216  return ( EXIT_FAILURE );
217 }
int main(int argc, char **argv)
MatrixEpetra< Real > matrix_Type
uint32_type ID
IDs.
Definition: LifeV.hpp:194
Real uTestFunctions(const Real &, const Real &x, const Real &y, const Real &z, const ID &i)
double Real
Generic real data.
Definition: LifeV.hpp:175
RegionMesh< LinearTetra > mesh_Type