LifeV
/timeAdvance.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 timeAdvance.cpp
28  \date 2010-02-10
29 
30  Author(s): F. Nobile <fabio.nobile@polimi.it>
31  M. Pozzoli <matteo1.pozzoli@mail.polimi.it>
32  C. Vergara <christian.vergara@polimi.it >
33  */
34 
35 /* ========================================================
36 
37 Simple Problem test with Dirichlet Boundary condition
38 
39 Solve the problem
40 
41  \frac{\partial u}{\partial t} - \Delta u = f
42 
43  u = u0 on the boundary
44 
45 linear_function.hpp:
46 
47  uexact = exp(-sin(Pi/2*t))*(x+y+z);
48 
49  f = (3 \pi^2 + 1 ) exp(-t) sin( \pi x) sin(\pi y) sin ( \pi z) on a cube
50 
51 nonlinear_function.hpp:
52 
53  uexact = exp(-sin(Pi/2*t))*cos(x *Pi)*cos(y*Pi)*cos(z*Pi);
54 
55  f = Pi2/4*( sin(Pi/2*t)+cos(Pi/2*t)*cos(Pi/2*t) )*exp(-sin(Pi/2*t))*cos(x *Pi)*cos(y*Pi)*cos(z*Pi);
56 */
57 
58 
59 // ===================================================
60 //! Includes
61 // ===================================================
62 
63 #include <Epetra_ConfigDefs.h>
64 #ifdef EPETRA_MPI
65 #include <mpi.h>
66 #include <Epetra_MpiComm.h>
67 #else
68 #include <Epetra_SerialComm.h>
69 #endif
70 
71 
72 #include <lifev/core/array/MapEpetra.hpp>
73 #include <lifev/core/mesh/MeshData.hpp>
74 #include <lifev/core/mesh/MeshPartitioner.hpp>
75 
76 #include <lifev/structure/solver/VenantKirchhoffViscoelasticSolver.hpp>
77 #include <lifev/structure/solver/VenantKirchhoffViscoelasticData.hpp>
78 
79 #include <lifev/core/fem/TimeData.hpp>
80 #include <lifev/core/fem/FESpace.hpp>
81 
82 #include <lifev/core/fem/TimeAdvance.hpp>
83 #include <lifev/core/fem/TimeAdvanceNewmark.hpp>
84 #include <lifev/core/fem/TimeAdvanceBDF.hpp>
85 
86 #include <lifev/core/filter/ExporterEnsight.hpp>
87 #include <lifev/core/filter/ExporterHDF5.hpp>
88 #include <lifev/core/filter/ExporterEmpty.hpp>
89 
90 #include <iostream>
91 
92 #include "timeAdvance.hpp"
93 
94 //choose of function.hpp
95 #include "linear_function.hpp"
96 //#include "nonlinear_function.hpp"
97 
98 // ===================================================
99 //! Namespaces & Define
100 // ===================================================
101 using namespace LifeV;
102 
103 const int TOP = 6;
104 const int BOTTOM = 5;
105 const int LEFT = 3;
106 const int RIGHT = 4;
107 const int FRONT = 2;
108 const int BACK = 1;
109 
110 
111 // ===================================================
112 //! Private members
113 // ===================================================
114 struct problem::Private
115 {
117  rho (1)
118  {}
119 
120  typedef std::function<Real ( Real const&, Real const&, Real const&, Real const&, ID const& ) > fct_type;
121  double rho;
123 
125 
127  {
128  fct_type f;
129  f = std::bind (&UZero, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3, std::placeholders::_4, std::placeholders::_5);
130  return f;
131  }
132 };
133 
134 
135 
136 
137 
138 // ===================================================
139 //! Constructors
140 // ===================================================
141 problem::problem ( int argc,
142  char** argv,
143  std::shared_ptr<Epetra_Comm> structComm) :
144  members ( new Private() )
145 {
146  GetPot command_line (argc, argv);
147  std::string data_file_name = command_line.follow ("data", 2, "-f", "--file");
148  GetPot dataFile ( data_file_name );
149  members->data_file_name = data_file_name;
150 
151  members->rho = dataFile ( "problem/physics/density", 1. );
152 
153  std::cout << "density = " << members->rho << std::endl;
154 
155  members->comm = structComm;
156  int ntasks = members->comm->NumProc();
157 
158  if (!members->comm->MyPID() )
159  {
160  std::cout << "My PID = " << members->comm->MyPID() << " out of " << ntasks << " running." << std::endl;
161  }
162 }
163 
164 // ===================================================
165 //! Methods
166 // ===================================================
167 void
168 problem::run()
169 {
170  typedef RegionMesh<LinearTetra> mesh_Type;
171  typedef VenantKirchhoffViscoelasticSolver< mesh_Type >::vector_type vector_type;
172  typedef std::shared_ptr<vector_type> vector_ptrtype;
173 
174  typedef std::shared_ptr< TimeAdvance< vector_type > > TimeAdvance_type;
175  typedef FESpace< mesh_Type, MapEpetra > FESpace_type;
176  typedef std::shared_ptr<FESpace_type> FESpace_ptrtype;
177 
178  bool verbose = (members->comm->MyPID() == 0);
179 
180  //
181  // VenantKirchhoffViscoelasticData
182  //
183 
184  GetPot dataFile ( members->data_file_name.c_str() );
185  std::shared_ptr<VenantKirchhoffViscoelasticData> dataProblem (new VenantKirchhoffViscoelasticData( ) );
186  dataProblem->setup (dataFile, "problem");
187 
188  MeshData meshData;
189  meshData.setup (dataFile, "problem/space_discretization");
190 
191  std::shared_ptr<mesh_Type > fullMeshPtr ( new mesh_Type ( members->comm ) );
192  readMesh (*fullMeshPtr, meshData);
193 
194  std::shared_ptr<mesh_Type > localMeshPtr;
195  {
196  MeshPartitioner< mesh_Type > meshPart ( fullMeshPtr, members->comm );
197  localMeshPtr = meshPart.meshPartition();
198  }
199 
200  //
201  // The Problem Solver
202  //
203 
204  // Scalar Solution Space:
205 
206  std::string Order = dataFile ( "problem/space_discretization/order", "P1");
207 
208  if ( Order.compare ("P1") == 0 )
209  {
210  if (verbose)
211  {
212  std::cout << " Space order : P1" << std::flush;
213  }
214 
215  }
216  else if ( Order.compare ("P2") == 0 )
217  {
218  if (verbose)
219  {
220  std::cout << " Space order : P2";
221  }
222  }
223  if (verbose)
224  {
225  std::cout << std::endl;
226  }
227 
228 
229  // finite element space of the solution
230 
231  std::string dOrder = dataFile ( "problem/space_discretization/order", "P1");
232 
233  FESpace_ptrtype feSpace ( new FESpace_type (localMeshPtr, dOrder, 1, members->comm) );
234 
235  // instantiation of the VenantKirchhoffViscoelasticSolver class
236 
237  VenantKirchhoffViscoelasticSolver< mesh_Type > problem;
238 
239  problem.setup (dataProblem, feSpace,
240  members->comm);
241 
242  problem.setDataFromGetPot (dataFile);
243 
244  // the boundary conditions
245 
246  BCFunctionBase uZero ( members->getUZero() );
247  BCFunctionBase uEx (uexact);
248 
249  BCHandler bcH;
250 
251  bcH.addBC ( "Top", TOP, Essential, Full, uEx, 1 );
252  bcH.addBC ( "Bottom", BOTTOM, Essential, Full, uEx, 1 );
253  bcH.addBC ( "Left", LEFT, Essential, Full, uEx, 1 );
254  bcH.addBC ( "Right", RIGHT, Essential, Full, uEx, 1 );
255  bcH.addBC ( "Front", FRONT, Essential, Full, uEx, 1 );
256  bcH.addBC ( "Back", BACK, Essential, Full, uEx, 1 );
257 
258  std::ofstream out_norm;
259  if (verbose)
260  {
261  out_norm.open ("norm.txt");
262  out_norm << " time "
263  << " L2_Error "
264  << " H1_Error "
265  << " L2_RelError "
266  << " H1_RelError \n";
267  out_norm.close();
268  }
269 
270  LifeChrono chrono;
271 
272  std::string TimeAdvanceMethod = dataFile ( "problem/time_discretization/method", "Newmark");
273 
274  TimeAdvance_type timeAdvance ( TimeAdvanceFactory::instance().createObject ( TimeAdvanceMethod ) );
275 
276  UInt OrderDev = 2;
277 
278  //! initialization of parameters of time Advance method:
279  if (TimeAdvanceMethod == "Newmark")
280  {
281  timeAdvance->setup ( dataProblem->dataTimeAdvance()->coefficientsNewmark() , OrderDev);
282  }
283 
284  if (TimeAdvanceMethod == "BDF")
285  {
286  timeAdvance->setup (dataProblem->dataTimeAdvance()->orderBDF() , OrderDev);
287  }
288 
289  timeAdvance->setTimeStep (dataProblem->dataTime()->timeStep() );
290  timeAdvance->showMe();
291 
292 
293  dataProblem->setup (dataFile, "problem");
294 
295  chrono.start();
296 
297  double xi = timeAdvance->coefficientSecondDerivative ( 0 ) / ( dataProblem->dataTime()->timeStep() * dataProblem->dataTime()->timeStep() );
298 
299  problem.buildSystem (xi);
300 
301 
302  MPI_Barrier (MPI_COMM_WORLD);
303 
304  if (verbose )
305  {
306  std::cout << "ok." << std::endl;
307  }
308 
309  // building some vectors:
310  MapEpetra uMap = problem.solution()->map();
311 
312  // computing the rhs
313  vector_type rhs ( uMap, Unique );
314 
315 
316  // postProcess
317  std::shared_ptr< Exporter<mesh_Type > > exporter;
318 
319  std::string const exporterType = dataFile ( "exporter/type", "ensight");
320 
321 #ifdef HAVE_HDF5
322  if (exporterType.compare ("hdf5") == 0)
323  {
324  exporter.reset ( new ExporterHDF5<mesh_Type > ( dataFile, "problem" ) );
325  }
326  else
327 #endif
328  {
329  if (exporterType.compare ("none") == 0)
330  {
331  exporter.reset ( new ExporterEmpty<mesh_Type > ( dataFile, localMeshPtr, "problem", members ->comm->MyPID() ) );
332  }
333  else
334  {
335  exporter.reset ( new ExporterEnsight<mesh_Type > ( dataFile, localMeshPtr, "problem", members->comm->MyPID() ) );
336  }
337  }
338  exporter->setPostDir ( "./" ); // This is a test to see if M_post_dir is working
339  exporter->setMeshProcId ( localMeshPtr, members->comm->MyPID() );
340 
341  vector_ptrtype U ( new vector_type (*problem.solution(), exporter->mapType() ) );
342  vector_ptrtype V ( new vector_type (*problem.solution(), exporter->mapType() ) );
343  vector_ptrtype W ( new vector_type (*problem.solution(), exporter->mapType() ) );
344  vector_ptrtype Exact ( new vector_type (*problem.solution(), exporter->mapType() ) );
345  vector_ptrtype vExact ( new vector_type (*problem.solution(), exporter->mapType() ) );
346  vector_ptrtype wExact ( new vector_type (*problem.solution(), exporter->mapType() ) );
347  vector_ptrtype RHS ( new vector_type (*problem.solution(), exporter->mapType() ) );
348  exporter->addVariable ( ExporterData<mesh_Type>::ScalarField, "displacement",
349  feSpace, U, UInt (0) );
350 
351  exporter->addVariable ( ExporterData<mesh_Type>::ScalarField, "velocity",
352  feSpace, V, UInt (0) );
353  exporter->addVariable ( ExporterData<mesh_Type>::ScalarField, "uexact",
354  feSpace, Exact, UInt (0) );
355 
356 
357  exporter->postProcess ( 0 );
358 
359  //evaluate disp and vel as interpolate the bcFunction d0 and v0
360  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( d0 ), *U, 0.0 );
361  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( v0 ), *V, 0.0 );
362  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( a0 ), *W, 0.0 );
363 
364  //evaluate disp and vel as interpolate the bcFunction d0 and v0 and w0
365 
366  std::vector<vector_ptrtype> uv0;
367 
368  Real dt = dataProblem->dataTime()->timeStep();
369  Real T = dataProblem->dataTime()->endTime();
370 
371  if (TimeAdvanceMethod == "Newmark")
372  {
373  uv0.push_back (U);
374  uv0.push_back (V);
375  uv0.push_back (W);
376  }
377  if (TimeAdvanceMethod == "BDF")
378  {
379  for ( UInt previousPass = 0; previousPass < dataProblem->dataTimeAdvance()->orderBDF() ; previousPass++)
380  {
381  Real previousTimeStep = -previousPass * dt;
382  feSpace->interpolate (static_cast<FESpace_type::function_Type> ( uexact ), *U, previousTimeStep );
383  uv0.push_back (U);
384  }
385  }
386 
387  timeAdvance->setInitialCondition (uv0);
388 
389  timeAdvance-> setTimeStep (dataProblem->dataTime()->timeStep() );
390 
391  timeAdvance->updateRHSContribution (dataProblem->dataTime()->timeStep() );
392 
393  vector_type uComputed (uMap, Repeated );
394  vector_type uExa (uMap, Repeated );
395  vector_type vExa (uMap, Repeated );
396  vector_type wExa (uMap, Repeated );
397 
398  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( uexact ), *Exact, 0 );
399  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( v0 ), *vExact, 0 );
400  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( a0 ), *wExact, 0 );
401 
402  *U = timeAdvance->solution();
403  *V = timeAdvance->firstDerivative();
404  *W = timeAdvance->secondDerivative();
405 
406 
407  exporter->postProcess ( 0 );
408 
409 
410 
411 
412  for (Real time = dt; time <= T; time += dt)
413  {
414  dataProblem->dataTime()->setTime (time);
415 
416  if (verbose)
417  {
418  std::cout << std::endl;
419  std::cout << " P - Now we are at time " << dataProblem->dataTime()->time() << " s." << std::endl;
420  }
421 
422  rhs *= 0;
423 
424  timeAdvance->updateRHSContribution ( dt );
425 
426  //evaluate rhs
427 
428  feSpace->l2ScalarProduct (source_in, rhs, time);
429  rhs += problem.matrMass() * timeAdvance->rhsContributionSecondDerivative();
430 
431  //update system
432  problem.updateRHS (rhs );
433 
434  //solver system
435  problem.iterate ( bcH ); // Computes the matrices and solves the system
436 
437  //update unknowns of timeAdvance
438  timeAdvance->shiftRight (*problem.solution() );
439 
440  //evaluate uexact solution
441  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( uexact ), *Exact, time );
442  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( v0 ), *vExact, time );
443  feSpace->interpolate ( static_cast<FESpace_type::function_Type> ( a0 ), *wExact, time );
444 
445  *U = timeAdvance->solution();
446  *V = timeAdvance->firstDerivative();
447  *W = timeAdvance->secondDerivative();
448 
449  //postProcess
450  exporter->postProcess ( time );
451 
452  // Error L2 and H1 Norms
453 
454  vector_type u (*problem.solution(), Repeated);
455 
456  AnalyticalSol uExact;
457 
458  Real H1_Error, H1_RelError, L2_Error, L2_RelError;
459 
460  L2_Error = feSpace->l2Error (uexact, u, time , &L2_RelError);
461  H1_Error = feSpace->h1Error (uExact, u, time , &H1_RelError);
462 
463  //save the norm
464  out_norm.open ("norm.txt", std::ios::app);
465  out_norm << time << " "
466  << L2_Error << " "
467  << H1_Error << " "
468  << L2_RelError << " "
469  << H1_RelError << "\n";
470  out_norm.close();
471 
472  MPI_Barrier (MPI_COMM_WORLD);
473  }
474  chrono.stop();
475  if (verbose)
476  {
477  std::cout << "Total iteration time " << chrono.diff() << " s." << std::endl;
478  }
479 }
FESpace - Short description here please!
Definition: FESpace.hpp:78
GetPot(const int argc_, char **argv_, const char *FieldSeparator=0x0)
Definition: GetPot.hpp:507
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
std::shared_ptr< std::vector< Int > > M_isOnProc
const int BOTTOM
MeshData - class for handling spatial discretization.
Definition: MeshData.hpp:72
std::function< Real(Real const &, Real const &, Real const &, Real const &, ID const &) > fct_type
double Real
Generic real data.
Definition: LifeV.hpp:175
std::shared_ptr< Epetra_Comm > comm
const int FRONT
fct_type getUZero()
Private members.
const int LEFT
const int BACK
const int RIGHT
std::string data_file_name
const int TOP
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191