LifeV
FSIMonolithicGI.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  *!
29  * @file
30  * @brief File containing the Monolithic Geometry--Implicit FSI Solver
31  *
32  * @date 18-09-2008
33  * @author Paolo Crosetto <paolo.crosetto@epfl.ch>
34  *
35  * @maintainer Paolo Crosetto <paolo.crosetto@epfl.ch>
36  */
37 
38 #ifndef _MONOLITHICGI_HPP
39 #define _MONOLITHICGI_HPP
40 
41 #include <lifev/fsi/solver/FSIMonolithic.hpp>
42 #include <lifev/core/array/MatrixBlockMonolithicEpetra.hpp>
43 #include <lifev/core/array/VectorBlockMonolithicEpetra.hpp>
44 
45 
46 namespace LifeV
47 {
48 #ifdef OBSOLETE
50 #endif
51 
52 typedef FactorySingleton< Factory< FSIOperator, std::string > > FSIFactory_Type;
53 
54 //! FSIMonolithicGI Geometry-Implicit solver
55 /*!
56  * @author Paolo Crosetto <paolo.crosetto@epfl.ch>
57  * Class handling the nonlinear monolithic solver for FSI problems. The (exact or inexact)
58  * Newton algorithm is used to solve the nonlinearity.
59  * The block structure of the jacobian matrix is
60  * \f$\left(\begin{array}{ccc} C&B&S\\ D&N&0\\ 0&E&H \end{array}\right)\f$
61  * where \f$N\f$ represents the solid block, \f$C\f$ the fluid block, \f$H\f$ is the harmonic extension block,
62  * while the extra
63  * diagonal blocks represent the coupling. The implementation of the stress continuity coupling condition
64  * is obtained by means of an augmented formulation.
65  * Different possible preconditioners are implemented.
66  *
67  * Important parameters to set properly in the data file:
68  * - useShapeDerivatives: if true the shape derivatives block is added to the Jacobian matrix;
69  * - domainVelImplicit: if true the domain velocity w in the convective term is considered an unknown (at the time n+1);
70  * - convectiveTermDer: false if the convective term is linearized (\f$u^{n+1}\nabla(u^n-w^n)\f$),
71  * otherwise it can be either true (if we use the Newton method to solve the convective term nonlinearity) or false
72  * (fixed-point method);
73  * - semiImplicit: if true only one iteration of the nonlinear solver is performed. Otherwise
74  * the nonlinear iterations continue up to the specified tolerance. Set it to true for the GCE;
75  * - method: can be either monolithicGE, monolithicGI if the geometry is treated respectively explicitly or implicitly,
76  * or exactJacobians, fixedPoint for partitioned strategies;
77  * - blockOper: specifies the matrix type to be used for the linear system: if AdditiveSchwarz, the matrix is the standard
78  * ine for GE; if AdditiveSchwarzRN the coupling blocks are of Robin type instead of Dirichlet and Neumann. The parameters
79  * for the Robin coupling are alphaf and alphas in the data file. NOTE: this method has currently been tested only for
80  * alphas=0.
81  * - DDBlockPrec: specifies the possible preconditioners to use. Can be: AdditiveSchwarz, MonolithicBlockComposedDN, MonolithicBlockComposedDN2,
82  * MonolithicBlockComposedNN, MonolithicBlockComposedDNND.
83  */
84 
85 class FSIMonolithicGI : public FSIMonolithic
86 {
87 public:
88 
89  typedef FSIMonolithic super_Type;
90  typedef Preconditioner prec_Type;
91  typedef std::shared_ptr< prec_Type > prec_type;
92 
93  //!@name Constructor and Destructor
94  //@{
95 
96  //! Empty Constructor
97  FSIMonolithicGI();
98 
99  //! Destructor
100  virtual ~FSIMonolithicGI() {}
101 
102  //@}
103 
104  //!@name Public Methods
105  //@{
106 
107  /**
108  Sets the parameters read from data file
109  */
110  void setup ( const GetPot& dataFile );
111 
112  //! initializes the fluid and mesh problems, creates the map of the global matrix
113  void setupFluidSolid ( UInt const fluxes );
114 
115  //! builds the constant part of the fluid-structure-mesh motion matrix
116  void buildSystem();
117 
118  /**
119  evaluates the residual b-Ax
120  \param res: output
121  \param _sol: fluid domain displacement solution
122  \param iter: current NonLinearRichardson (block Gauss Seidel for the tangent system) iteration
123  */
124  void evalResidual ( vector_Type& res, const vector_Type& sol, const UInt iter );
125 
126  //!Apply the boundary conditions to each block composing the monolithic problem
127  /**
128  Sets the vectors of: boundary conditions, FESpaces, couplings, offsets, and sets the blocks in the composed operator
129  which constitutes the monolithic problem. Then calls the applyBoundaryConditions of the MonolithicBlockMatrix operator, passing
130  also the right hand side.
131  */
132  void applyBoundaryConditions();
133 
134  void updateSolution ( const vector_Type& solution )
135  {
136  super_Type::updateSolution ( solution );
137 
138  //The size of the vectors for the ALE is = dimension of the ALE problem
139  //To do the shift right we first need to extract the fluid displacement
140  //And then push it into the ALE timeAdvance class.
141  vectorPtr_Type displacementToSave ( new vector_Type (M_mmFESpace->map() ) );
142  UInt offset ( M_solidAndFluidDim + nDimensions * M_interface );
143  displacementToSave->subset (solution, offset);
144 
145  //This updateRHSFirstDerivative has to be done before the shiftRight
146  //In fact it updates the right hand side of the velocity using the
147  //previous times. The method velocity() uses it and then, the compuation
148  //of the velocity is done using the current time and the previous times.
149  //M_ALETimeAdvance->updateRHSFirstDerivative ( M_data->dataFluid()->dataTime()->timeStep() );
150  M_ALETimeAdvance->shiftRight ( *displacementToSave );
151  }
152 
153  //! Set vectors for restart
154  /*!
155  * Set vectors for restart
156  */
157  void setALEVectorInStencil (const vectorPtr_Type& fluidDisp,
158  const UInt iter,
159  const bool lastVector);
160 
161  //!@name Get Methods
162  //@{
163 
164  //! get the current solution vector.
165  LIFEV_DEPRECATED ( const vector_Type& solution() const )
166  {
167  if ( M_epetraWorldComm->MyPID() == 0 )
168  {
169  std::cerr << std::endl << "Warning: FSIMonolithic::solution() is deprecated!" << std::endl
170  << " You should not access the solution inside FSIOperator or FSIMonolithic!" << std::endl;
171  }
172 
173  return M_fluidTimeAdvance->singleElement (0);
174  }
175 
176  //! getter for the map of fluid-structure-interface (without the mesh motion)
177  const MapEpetra& mapWithoutMesh() const
178  {
179  return *M_mapWithoutMesh;
180  }
181 
182  //! getter for the global matrix of the system
183  const matrixPtr_Type matrixPtr() const
184  {
185  return M_monolithicMatrix->matrix();
186  }
187 
188  static bool S_register;
189 
190  //@}
191 
192 protected:
193 
194  //!@name Protected Methods
195  //@{
196 
197  //! set the block preconditioner
198  void setupBlockPrec();
199 
200  //@}
201 
202 private:
203 
204  //! @name Private Methods
205  //@{
206 
207  //! Factory method for the system matrix, of type MonolithicBlockBase
208  void createOperator ( std::string& operType )
209  {
210  M_monolithicMatrix.reset ( MonolithicBlockMatrix::Factory_Type::instance().createObject ( operType ) );
211  M_monolithicMatrix.reset ( MonolithicBlockMatrix::Factory_Type::instance().createObject ( operType ) );
212  }
213 
214  /**
215  calculates the terms due to the shape derivatives given the mesh increment deltaDisp. The shape derivative block is assembled in a matrix
216  (not in a right hand side representing the matrix-vector multiplication)
217  \param sdMatrix: output. Shape derivatives block to be summed to the Jacobian matrix.
218  */
219  void shapeDerivatives ( FSIOperator::fluid_Type::matrixPtr_Type sdMatrix );
220 
221  //! assembles the mesh motion matrix.
222  /*!In Particular it diagonalize the part of the matrix corresponding to the
223  Dirichlet condition expressing the coupling
224  \param iter: current iteration: used as flag to distinguish the first nonlinear iteration from the others
225  */
226  void assembleMeshBlock ( UInt iter );
227 
228  //@}
229 
230 
231  //!@name Private Members
232  //@{
233 
234  std::shared_ptr<MapEpetra> M_mapWithoutMesh;
235  //This vector is used in the shapeDerivatives method since a
236  //copy of the solution at the current iteration k is necessary
237  vectorPtr_Type M_uk;
238  UInt M_interface;
239  matrixPtr_Type M_meshBlock;
240  FSIOperator::fluid_Type::matrixPtr_Type M_shapeDerivativesBlock;
241  matrixPtr_Type M_solidDerBlock;
242  //std::vector<fluidBchandlerPtr_Type> M_BChsLin;
243  //@}
244 
245  //! Factory method
246  static FSIOperator* instantiate()
247  {
248  return new FSIMonolithicGI();
249  }
250 
251 };
252 
253 //! Factory create function
254 inline FSIMonolithic* createFSIMonolithicGI()
255 {
256  return new FSIMonolithicGI();
257 }
258 
259 }
260 #endif
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
const UInt nDimensions(NDIM)
#define LIFEV_DEPRECATED(func)
Definition: LifeV.hpp:117
uint32_type UInt
generic unsigned integer (used mainly for addressing)
Definition: LifeV.hpp:191