LifeV
FSIData.hpp
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 FSIData - File containing a data container for FSI problems
30 
31  @author Cristiano Malossi <cristiano.malossi@epfl.ch>
32  @author Gilles fourestey <gilles.fourestey@epfl.ch>
33  @date 10-06-2010
34 
35  @maintainer Simone Deparis <simone.deparis@epfl.ch>
36  */
37 
38 
39 #ifndef FSIDATA_H
40 #define FSIDATA_H
41 
42 #include <lifev/navier_stokes/solver/OseenData.hpp>
43 #include <lifev/structure/solver/StructuralConstitutiveLawData.hpp>
44 
45 #include <array>
46 
47 namespace LifeV
48 {
49 
50 //! FSIData - Data container for FSI problems
51 /*!
52  * @author Cristiano Malossi
53  */
54 
55 class FSIData
56 {
57 public:
58 
59  //! @name Type definitions
60  //@{
61 
62  typedef OseenData dataFluid_Type;
63  typedef std::shared_ptr< dataFluid_Type > dataFluidPtr_Type;
64 
65  typedef StructuralConstitutiveLawData dataSolid_Type;
66  typedef std::shared_ptr< dataSolid_Type > dataSolidPtr_Type;
67 
68  typedef dataFluid_Type::time_Type time_Type;
69  typedef dataFluid_Type::timePtr_Type timePtr_Type;
70 
71  typedef dataFluid_Type::timeAdvance_Type timeAdvance_Type;
72  typedef dataFluid_Type::timeAdvancePtr_Type timeAdvancePtr_Type;
73 
74  //@}
75 
76 
77  //! @name Constructors & Destructor
78  //@{
79 
80  //! Empty Constructor
81  explicit FSIData();
82 
83  //! Copy constructor
84  /*!
85  * @param FSIData - FSIData
86  */
87  explicit FSIData ( const FSIData& FSIData );
88 
89  //! Destructor
90  virtual ~FSIData() {}
91 
92  //@}
93 
94 
95  //! @name Operators
96  //@{
97 
98  //! Operator=
99  /*!
100  * @param FSIData - FSIData
101  */
102  FSIData& operator= ( const FSIData& FSIData );
103 
104  //@}
105 
106 
107  //! @name Methods
108  //@{
109 
110  //! Read the dataFile and set all the quantities
111  /*!
112  * @param dataFile - data file
113  */
114  void setup ( const GetPot& dataFile, const std::string& section = "problem" );
115 
116  //! Display the values
117  void showMe ( std::ostream& output = std::cout );
118 
119  bool isMonolithic();
120 
121  //@}
122 
123 
124  //! @name Set methods
125  //@{
126 
127  //! Set data fluid container
128  /*!
129  * @param dataFluid shared_ptr to dataFluid container
130  */
131  void setDataFluid ( const dataFluidPtr_Type& dataFluid )
132  {
133  M_dataFluid = dataFluid;
134  }
135 
136  //! Set data solid container
137  /*!
138  * @param dataFluid shared_ptr to dataSolid container
139  */
140  void setDataSolid ( const dataSolidPtr_Type& dataSolid )
141  {
142  M_dataSolid = dataSolid;
143  }
144 
145 
146  //! Set TimeData ALE container
147  /*!
148  * @param timeDataALE shared_ptr to timeDataALE container
149  */
150  void setTimeDataALE ( const timePtr_Type& timeALE )
151  {
152  M_timeALE = timeALE;
153  }
154 
155  //! Set TimeData ALE container
156  /*!
157  * @param timeDataALE shared_ptr to timeDataALE container
158  */
159  void setTimeAdvanceDataALE ( const timeAdvancePtr_Type& timeAdvanceALE )
160  {
161  M_timeAdvanceALE = timeAdvanceALE;
162  }
163 
164  //@}
165 
166 
167  //! @name Get methods
168  //@{
169 
170  //! Get data fluid container
171  /*!
172  * @return shared_ptr to dataFluid container
173  */
174  const dataFluidPtr_Type& dataFluid() const
175  {
176  return M_dataFluid;
177  }
178 
179  //! Get data solid container
180  /*!
181  * @return shared_ptr to dataSolid container
182  */
183  const dataSolidPtr_Type& dataSolid() const
184  {
185  return M_dataSolid;
186  }
187 
188  //! Get data time ALE
189  /*!
190  * @return shared_ptr to ALE dataTime container
191  */
192  const timePtr_Type& timeDataALE() const
193  {
194  return M_timeALE;
195  }
196 
197  //! Get data time ALE
198  /*!
199  * @return shared_ptr to ALE dataTimeAdvance container
200  */
201  const timeAdvancePtr_Type& timeAdvanceDataALE() const
202  {
203  return M_timeAdvanceALE;
204  }
205 
206  //! Get maximum number of subiterations
207  /*!
208  * @return maximum number of subiterations
209  */
210  const UInt& maxSubIterationNumber() const
211  {
212  return M_maxSubIterationNumber;
213  }
214 
215  //! Get absolute tolerance
216  /*!
217  * @return absolute tolerance
218  */
219  const Real& absoluteTolerance() const
220  {
221  return M_absoluteTolerance;
222  }
223 
224  //! Get relative tolerance
225  /*!
226  * @return relative tolerance
227  */
228  const Real& relativeTolerance() const
229  {
230  return M_relativeTolerance;
231  }
232 
233  //! Get error tolerance
234  /*!
235  * @return error tolerance
236  */
237  const Real& errorTolerance() const
238  {
239  return M_errorTolerance;
240  }
241 
242  //! Get NonLinearLineSearch
243  /*!
244  * @return NonLinearLineSearch
245  */
246  const Int& NonLinearLineSearch() const
247  {
248  return M_NonLinearLineSearch;
249  }
250 
251  //! Get method type
252  /*!
253  * @return method type
254  */
255  const std::string& method() const
256  {
257  return M_method;
258  }
259 
260  //! Get algorithm type
261  /*!
262  * @return algorithm type
263  */
264  const std::string& algorithm() const
265  {
266  return M_algorithm;
267  }
268 
269  //! Get default omega for Aitken iterations
270  /*!
271  * @return default omega for Aitken iterations
272  */
273  const Real& defaultOmega() const
274  {
275  return M_defaultOmega;
276  }
277 
278  //! Get the range of omega for Aitken iterations
279  /*!
280  * @return range of omega for Aitken iterations
281  */
282  const std::array< Real, 2 >& OmegaRange() const
283  {
284  return M_rangeOmega;
285  }
286 
287  //! Get update every
288  /*!
289  * If M_updateEvery == 1, normal fixedPoint algorithm
290  * If M_updateEvery > 1, recompute computational domain every M_updateEvery iterations (transpiration)
291  * If M_updateEvery <= 0, recompute computational domain and matrices only at first subiteration (semi-implicit)
292  *
293  * @return updateEvery value
294  */
295  const Int& updateEvery() const
296  {
297  return M_updateEvery;
298  }
299 
300  //! Get the fluid Interface Flag
301  /*!
302  * @return Flag of the interface on the fluid boundary side
303  */
304  const Int& fluidInterfaceFlag() const
305  {
306  return M_fluidInterfaceFlag;
307  }
308 
309  //! Get the structure Interface Flag
310  /*!
311  * @return Flag of the interface on the structure boundary side
312  */
313  const Int& structureInterfaceFlag() const
314  {
315  return M_structureInterfaceFlag;
316  }
317 
318  //! Get the fluid Interface Flag (for Vertices)
319  /*!
320  * @return Flag of the vertex on the interface on the fluid boundary side
321  */
322  Int const* fluidInterfaceVertexFlag() const
323  {
324  return M_fluidInterfaceVertexFlag.get();
325  }
326 
327  //! Get the fluid Interface Flag (for Vertices)
328  /*!
329  * @return Flag of the vertex on the interface on the structure boundary side
330  */
331  Int const* structureInterfaceVertexFlag() const
332  {
333  return M_structureInterfaceVertexFlag.get();
334  }
335 
336  //! Get the tolerance for the Interface identification
337  /*!
338  * @return the tolerance for the Interface identification
339  */
340  const Real& interfaceTolerance() const
341  {
342  return M_interfaceTolerance;
343  }
344 
345  //! Get the timestep to restart the simulation
346  /*!
347  * @return the timestep used in the previous simulation from which we want to restart, used for the initialization
348  of the time discretization
349  */
350  inline Real restartTimeStep() const
351  {
352  return M_restartTimeStep;
353  }
354 
355 
356  //@}
357 
358 private:
359 
360  dataFluidPtr_Type M_dataFluid;
361  dataSolidPtr_Type M_dataSolid;
362  timePtr_Type M_timeALE;
363  timeAdvancePtr_Type M_timeAdvanceALE;
364 
365  // Problem - Non Linear Richardson parameters
366  UInt M_maxSubIterationNumber;
367  Real M_absoluteTolerance;
368  Real M_relativeTolerance;
369  Real M_errorTolerance;
370  Int M_NonLinearLineSearch;
371 
372  // Problem - Methods
373  std::string M_method;
374  std::string M_algorithm;
375 
376  // Problem - FixedPoint / EJ
377  Real M_defaultOmega;
378  std::array< Real, 2 > M_rangeOmega;
379  Int M_updateEvery;
380 
381  // Interface
382  Int M_fluidInterfaceFlag;
383  Int M_structureInterfaceFlag;
384 
385  std::unique_ptr<Int const> M_fluidInterfaceVertexFlag;
386  std::unique_ptr<Int const> M_structureInterfaceVertexFlag;
387 
388  Real M_interfaceTolerance;
389 
390  Real M_restartTimeStep;
391 };
392 
393 } // end namespace LifeV
394 
395 #endif // end FSIDATA_H
std::shared_ptr< timeAdvance_Type > timeAdvancePtr_Type
Definition: OseenData.hpp:110
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
TimeAdvanceData timeAdvance_Type
Definition: OseenData.hpp:109
int32_type Int
Generic integer data.
Definition: LifeV.hpp:188
std::shared_ptr< time_Type > timePtr_Type
Definition: OseenData.hpp:107
double Real
Generic real data.
Definition: LifeV.hpp:175
TimeData time_Type
Definition: OseenData.hpp:106
OseenData - LifeV Base class which holds usual data for the NavierStokes equations solvers...
Definition: OseenData.hpp:98
DataElasticStructure - Data container for solid problems with elastic structure.
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191