LifeV
lifev/fsi_blocks/testsuite/fsi_tube/main.cpp
Go to the documentation of this file.
1 /* -*- mode: c++ -*-
2 
3  This file is part of the LifeV library.
4  Copyright (C) 2010 EPFL
5 
6  This program is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2.1 of the License, or
9  (at your option) any later version.
10 
11  This program is distributed in the hope that it will be useful, but
12  WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14  General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with this program; if not, write to the Free Software
18  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
19  USA
20  */
21 /**
22  \file main.cpp
23  \author Davide Forti <davide.forti@epfl.ch>
24  \date 2014-02-06
25  */
26 
27 
28 #include <Epetra_ConfigDefs.h>
29 #ifdef EPETRA_MPI
30 #include <mpi.h>
31 #include <Epetra_MpiComm.h>
32 #else
33 #include <Epetra_SerialComm.h>
34 #endif
35 
36 #include <lifev/core/LifeV.hpp>
37 
38 #include <lifev/core/mesh/MeshData.hpp>
39 #include <lifev/core/mesh/MeshPartitioner.hpp>
40 
41 #include <lifev/core/array/VectorEpetra.hpp>
42 
43 #include <lifev/fsi_blocks/solver/FSIHandler.hpp>
44 
46 
47 using namespace LifeV;
48 
49 int
50 main ( int argc, char** argv )
51 {
52  bool verbose (false);
53 #ifdef HAVE_MPI
54  MPI_Init (&argc, &argv);
55  std::shared_ptr<Epetra_Comm> Comm ( new Epetra_MpiComm (MPI_COMM_WORLD) );
56  if ( Comm->MyPID() == 0 )
57  {
58  verbose = true;
59  }
60 #else
61  std::shared_ptr<Epetra_Comm> Comm( new Epetra_SerialComm () );
62  verbose = true;
63 #endif
64 
65  typedef RegionMesh<LinearTetra> mesh_Type;
66  typedef VectorEpetra vector_Type;
67  typedef std::shared_ptr<vector_Type> vectorPtr_Type;
68 
69  Real normTwo_sol;
70 
71  {
72  // ---------------------------//
73  // Reading the input datafile //
74  // ---------------------------//
75 
76  const std::string defaultDataName = "data";
77  GetPot command_line (argc, argv);
78  std::string data_file_name = command_line.follow (defaultDataName.c_str(), 2, "-f", "--file");
79  GetPot dataFile( data_file_name );
80 
81  // --------------------------------------------------//
82  // Initializing the handler and setting the datafile //
83  // --------------------------------------------------//
84 
85  FSIHandler fsi ( Comm );
86 
87  fsi.setDatafile ( dataFile );
88 
89  // check about the use of meshes: using partitioned meshes or do we have to partition them online?
90 
91  bool usePartitionedMeshes = dataFile("offlinePartioner/readPartitionedMeshes", false);
92 
93  if ( usePartitionedMeshes )
94  {
95  if ( verbose )
96  std::cout << "\n[PreProcessing] - Using partitioned meshes\n";
97 
98  int numParts = dataFile("offlinePartioner/parts", 2);
99 
100  if( numParts != Comm->NumProc())
101  {
102  if ( verbose )
103  std::cout << "\n[PreProcessing] - Using the wrong number of processes\n";
104 
105  return EXIT_FAILURE;
106  }
107 
108  // -------------------------------------------------------//
109  // Loading the partitioned fluid and the structure meshes //
110  // -------------------------------------------------------//
111 
112  fsi.readPartitionedMeshes ( );
113  }
114  else
115  {
116  if ( verbose )
117  std::cout << "\n[PreProcessing] - The meshes will be partioned during this simulation\n\n";
118 
119  // -------------------------------------------//
120  // Loading the fluid and the structure meshes //
121  // -------------------------------------------//
122 
123  fsi.readMeshes ( );
124 
125  // ------------------------------------------------//
126  // Partitioning the fluid and the structure meshes //
127  // ------------------------------------------------//
128 
129  fsi.partitionMeshes ( );
130  }
131 
132 
133  // --------------------------------//
134  // Getting the boundary conditions //
135  // --------------------------------//
136 
137  std::shared_ptr<BCHandler> interfaceFluidBC ( new BCHandler (*BCh_interfaceFluid ( ) ) );
138  fsi.setFluidInterfaceBoundaryConditions(interfaceFluidBC);
139 
140  std::shared_ptr<BCHandler> fluidBC ( new BCHandler (*BCh_fluid () ) );
141  std::shared_ptr<BCHandler> structureBC ( new BCHandler (*BCh_structure () ) );
142  std::shared_ptr<BCHandler> aleBC ( new BCHandler (*BCh_ale () ) );
143  fsi.setBoundaryConditions( fluidBC, structureBC, aleBC );
144 
145  std::string preconditioner = dataFile("fluid/preconditionerType","none");
146 
147  // ---------------------------------------------------------------------------------------------//
148  // Reading the physical informations for the fluid and the structure and initialize the solvers //
149  // ---------------------------------------------------------------------------------------------//
150 
151  fsi.setup();
152 
153  // ---------------------//
154  // Create inteface maps //
155  // ---------------------//
156 
157  fsi.buildInterfaceMaps ( );
158 
159  // -------------------------------//
160  // Create blocks for the coupling //
161  // -------------------------------//
162 
163  if ( !dataFile( "interface/nonconforming", false ) )
164  fsi.assembleCoupling ( );
165 
166  // ----------//
167  // Time Loop //
168  // ----------//
169 
170  fsi.solveFSIproblem ( );
171 
172  normTwo_sol = fsi.getFSIsolution()->norm2();
173 
174  }
175 
176 #ifdef HAVE_MPI
177  if (verbose)
178  {
179  std::cout << "\nMPI Finalization\n" << std::endl;
180  }
181  MPI_Finalize();
182 #endif
183 
184  if ( std::abs(normTwo_sol - 162041.4417 ) < 1.0e-3 )
185  {
186  return ( EXIT_SUCCESS );
187  }
188  else
189  {
190  return ( EXIT_FAILURE );
191  }
192 
193 }
GetPot(const int argc_, char **argv_, const char *FieldSeparator=0x0)
Definition: GetPot.hpp:507
double Real
Generic real data.
Definition: LifeV.hpp:175
int operator()(const char *VarName, int Default) const
Definition: GetPot.hpp:2021
std::shared_ptr< ExporterHDF5< mesh_Type > > M_importerStructure
Definition: FSIHandler.hpp:562
bool operator()(const char *VarName, bool Default) const
Definition: GetPot.hpp:2008
int main(int argc, char **argv)