26 #include <life/lifecore/LifeV.hpp> 27 #include <life/lifecore/LifeChrono.hpp> 28 #include <life/lifesolver/OseenSolver.hpp> 29 #include <life/lifemesh/MeshData.hpp> 30 #include <life/lifesolver/OseenData.hpp> 31 #include <life/lifearray/MapEpetra.hpp> 32 #include <life/lifefem/FESpace.hpp> 33 #include <life/lifefem/TimeAdvanceBDFNavierStokes.hpp> 34 #include <life/lifemesh/MeshPartitioner.hpp> 35 #include <life/lifefilters/ExporterEnsight.hpp> 36 #include <life/lifefilters/ExporterHDF5.hpp> 37 #include <life/lifefilters/ExporterEmpty.hpp> 41 #include "Epetra_config.h" 43 #include "Epetra_MpiComm.h" 45 #include "Epetra_SerialComm.h" 84 int main (
int argc,
char** argv)
91 MPI_Init (&argc, &argv);
94 std::shared_ptr<Epetra_Comm> comm;
96 comm.reset (
new Epetra_MpiComm ( MPI_COMM_WORLD ) );
98 MPI_Comm_size (MPI_COMM_WORLD, &nproc);
100 comm.reset (
new Epetra_SerialComm() );
103 bool verbose (
false);
104 if (comm->MyPID() == 0)
108 <<
" +-----------------------------------------------+" << std::endl
109 <<
" | Cavity example for LifeV |" << std::endl
110 <<
" +-----------------------------------------------+" << std::endl
112 <<
" +-----------------------------------------------+" << std::endl
113 <<
" | Author: Gwenol Grandperrin |" << std::endl
114 <<
" | Date: October 25, 2010 |" << std::endl
115 <<
" +-----------------------------------------------+" << std::endl
118 std::cout <<
"[Initilization of MPI]" << std::endl;
120 std::cout <<
"Using MPI (" << nproc <<
" proc.)" << std::endl;
122 std::cout <<
"Using serial version" << std::endl;
127 LifeV::Chrono globalChrono;
128 LifeV::Chrono initChrono;
129 LifeV::Chrono iterChrono;
130 globalChrono.start();
138 std::cout << std::endl <<
"[Loading the data]" << std::endl;
140 GetPot command_line (argc, argv);
141 const std::string dataFileName = command_line.follow (
"data", 2,
"-f",
"--file");
142 GetPot dataFile (dataFileName);
149 std::cout << std::endl <<
"[Loading the mesh]" << std::endl;
151 LifeV::MeshData meshData;
152 meshData.setup (dataFile,
"fluid/space_discretization");
155 std::cout <<
"Mesh file: " << meshData.meshDir() << meshData.meshFile() << std::endl;
157 std::shared_ptr< LifeV::RegionMesh<LifeV::LinearTetra> > fullMeshPtr (
new LifeV::RegionMesh<LifeV::LinearTetra>);
158 LifeV::readMesh (*fullMeshPtr, meshData);
160 LifeV::partitionMesh< LifeV::RegionMesh<LifeV::LinearTetra> > meshPart (fullMeshPtr, comm);
167 std::cout << std::endl <<
"[Creating the FE spaces]" << std::endl;
169 std::string uOrder = dataFile (
"fluid/space_discretization/vel_order",
"P2");
170 std::string pOrder = dataFile (
"fluid/space_discretization/press_order",
"P1");
171 if (verbose) std::cout <<
"FE for the velocity: " << uOrder << std::endl
172 <<
"FE for the pressure: " << pOrder << std::endl;
176 std::cout <<
"Building the velocity FE space... " << std::flush;
178 LifeV::FESpace< LifeV::RegionMesh<LifeV::LinearTetra>, LifeV::MapEpetra > uFESpace (meshPart, uOrder, 3, comm);
181 std::cout <<
"ok." << std::endl;
186 std::cout <<
"Building the pressure FE space... " << std::flush;
188 LifeV::FESpace< LifeV::RegionMesh<LifeV::LinearTetra>, LifeV::MapEpetra > pFESpace (meshPart, pOrder, 1, comm);
191 std::cout <<
"ok." << std::endl;
195 LifeV::UInt totalVelDof = uFESpace.map().getMap (LifeV::Unique)->NumGlobalElements();
196 LifeV::UInt totalPressDof = pFESpace.map().getMap (LifeV::Unique)->NumGlobalElements();
200 std::cout <<
"Total Velocity Dof: " << totalVelDof << std::endl;
204 std::cout <<
"Total Pressure Dof: " << totalPressDof << std::endl;
212 std::cout << std::endl <<
"[Boundary conditions]" << std::endl;
215 LifeV::BCFunctionBase uZero (fZero);
216 LifeV::BCFunctionBase uLid (lidBC);
218 LifeV::BCHandler::BCHints hint;
219 hint = LifeV::BCHandler::HINT_BC_NONE;
221 std::vector<LifeV::ID> xComp (1);
224 LifeV::BCHandler bcH ( 0, hint );
226 bcH.addBC (
"Top" , TOP , LifeV::Essential, LifeV::Full , uLid , 3 );
227 bcH.addBC (
"Left" , LEFT , LifeV::Essential, LifeV::Full , uZero, 3 );
228 bcH.addBC (
"Front" , FRONT , LifeV::Essential, LifeV::Component, uZero, xComp );
229 bcH.addBC (
"Right" , RIGHT , LifeV::Essential, LifeV::Full , uZero, 3 );
230 bcH.addBC (
"Back" , BACK , LifeV::Essential, LifeV::Component, uZero, xComp );
231 bcH.addBC (
"Bottom", BOTTOM, LifeV::Essential, LifeV::Full , uZero, 3 );
234 std::vector<LifeV::bcName_Type> fluxVector = bcH.getBCWithType ( LifeV::Flux );
235 LifeV::UInt numLM =
static_cast<LifeV::UInt> ( fluxVector.size() );
237 LifeV::UInt offset = uFESpace.map().getMap (LifeV::Unique)->NumGlobalElements()
238 + pFESpace.map().getMap (LifeV::Unique)->NumGlobalElements();
240 for ( LifeV::UInt i = 0; i < numLM; ++i )
242 bcH.setOffset ( fluxVector[i], offset + i );
250 std::cout << std::endl <<
"[Creating the problem]" << std::endl;
252 std::shared_ptr<LifeV::OseenData> oseenData (
new LifeV::OseenData() );
253 oseenData->setup ( dataFile );
257 std::cout <<
"Time discretization order " << oseenData->dataTime()->orderBDF() << std::endl;
261 LifeV::OseenSolver< LifeV::RegionMesh<LifeV::LinearTetra> > fluid (oseenData,
267 fluid.setUp (dataFile);
273 LifeV::MapEpetra fullMap (fluid.getMap() );
276 MPI_Barrier (MPI_COMM_WORLD);
283 std::cout << std::endl <<
"[Initialization of the simulation]" << std::endl;
285 LifeV::Real dt = oseenData->dataTime()->timeStep();
286 LifeV::Real t0 = oseenData->dataTime()->initialTime();
287 LifeV::Real tFinal = oseenData->dataTime()->endTime();
290 LifeV::TimeAdvanceBDFNavierStokes<vector_type> bdf (oseenData->dataTime()->orderBDF() );
293 t0 -= dt * bdf.bdfVelocity().order();
295 vector_type beta ( fullMap );
296 vector_type rhs ( fullMap );
298 MPI_Barrier (MPI_COMM_WORLD);
301 oseenData->dataTime()->setTime (t0);
305 fluid.updateSystem (0.0, beta, rhs);
307 bdf.bdfVelocity().setInitialCondition ( *fluid.solution() );
309 LifeV::Real time = t0 + dt;
310 for ( ; time <= oseenData->dataTime()->initialTime() + dt / 2.; time += dt)
312 oseenData->dataTime()->setTime (time);
314 fluid.updateSystem (0.0, beta, rhs);
316 bdf.bdfVelocity().shiftRight ( *fluid.solution() );
323 std::shared_ptr< LifeV::Exporter<LifeV::RegionMesh<LifeV::LinearTetra> > > exporter;
325 vector_ptrtype velAndPressure;
327 std::string
const exporterType = dataFile (
"exporter/type",
"ensight");
330 if (exporterType.compare (
"hdf5") == 0)
332 exporter.reset (
new LifeV::Hdf5exporter<LifeV::RegionMesh<LifeV::LinearTetra> > ( dataFile,
"cavity_example" ) );
333 exporter->setDirectory (
"./" );
334 exporter->setMeshProcId ( meshPart.meshPartition(), comm->MyPID() );
339 if (exporterType.compare (
"none") == 0)
341 exporter.reset (
new LifeV::NoExport<LifeV::RegionMesh<LifeV::LinearTetra> > ( dataFile, meshPart.meshPartition(),
"cavity_example", comm->MyPID() ) );
345 exporter.reset (
new LifeV::Ensight<LifeV::RegionMesh<LifeV::LinearTetra> > ( dataFile, meshPart.meshPartition(),
"cavity_example", comm->MyPID() ) );
349 velAndPressure.reset (
new vector_type (*fluid.solution(), exporter->mapType() ) );
351 exporter->addVariable ( LifeV::ExporterData::Vector,
"velocity", velAndPressure,
352 LifeV::UInt (0), uFESpace.dof().numTotalDof() );
354 exporter->addVariable ( LifeV::ExporterData::Scalar,
"pressure", velAndPressure,
355 LifeV::UInt (3 * uFESpace.dof().numTotalDof() ),
356 LifeV::UInt (pFESpace.dof().numTotalDof() ) );
357 exporter->postProcess ( 0 );
362 std::cout <<
"Initialization time: " << initChrono.diff() <<
" s." << std::endl;
370 std::cout << std::endl <<
"[Solving the problem]" << std::endl;
374 for ( ; time <= tFinal + dt / 2.; time += dt, iter++)
377 oseenData->dataTime()->setTime (time);
381 std::cout <<
"[t = " << oseenData->dataTime()->time() <<
" s.]" << std::endl;
386 double alpha = bdf.bdfVelocity().coefficientFirstDerivative ( 0 ) / oseenData->dataTime()->timeStep();
388 beta = bdf.bdfVelocity().extrapolation();
390 bdf.bdfVelocity().updateRHSContribution ( oseenData->dataTime()->timeStep() );
391 rhs = fluid.matrMass() * bdf.bdfVelocity().rhsContributionFirstDerivative();
392 fluid.getDisplayer().leaderPrint (
"alpha ", alpha);
393 fluid.getDisplayer().leaderPrint (
"\n");
394 fluid.getDisplayer().leaderPrint (
"norm beta ", beta.Norm2() );
395 fluid.getDisplayer().leaderPrint (
"\n");
396 fluid.getDisplayer().leaderPrint (
"norm rhs ", rhs.Norm2() );
397 fluid.getDisplayer().leaderPrint (
"\n");
399 fluid.updateSystem ( alpha, beta, rhs );
400 fluid.iterate ( bcH );
402 bdf.bdfVelocity().shiftRight ( *fluid.solution() );
405 vector_type vel (uFESpace.map(), LifeV::Repeated);
406 vector_type press (pFESpace.map(), LifeV::Repeated);
407 vector_type velpressure ( *fluid.solution(), LifeV::Repeated );
409 velpressure = *fluid.solution();
410 vel.subset (velpressure);
411 press.subset (velpressure, uFESpace.dim() *uFESpace.fieldDim() );
414 bool verbose = (comm->MyPID() == 0);
418 *velAndPressure = *fluid.solution();
419 exporter->postProcess ( time );
422 MPI_Barrier (MPI_COMM_WORLD);
427 std::cout <<
"Iteration time: " << iterChrono.diff() <<
" s." << std::endl << std::endl;
434 std::cout <<
"Total simulation time: " << globalChrono.diff() <<
" s." << std::endl;
LifeV::Real fZero(const LifeV::Real &, const LifeV::Real &, const LifeV::Real &, const LifeV::Real &, const LifeV::ID &)
fluid_type::vector_type vector_type
int main(int argc, char **argv)
LifeV::Real lidBC(const LifeV::Real &t, const LifeV::Real &, const LifeV::Real &, const LifeV::Real &, const LifeV::ID &i)
LifeV::OseenSolver< mesh_type > fluid_type
LifeV::RegionMesh< LifeV::LinearTetra > mesh_type
std::shared_ptr< vector_type > vector_ptrtype