LifeV
MonolithicBlockComposedDN.hpp
Go to the documentation of this file.
1 /* -*- mode: c++ -*- */
2 //@HEADER
3 /*
4 *******************************************************************************
5 
6  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
7  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
8 
9  This file is part of LifeV.
10 
11  LifeV is free software; you can redistribute it and/or modify
12  it under the terms of the GNU Lesser General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  LifeV is distributed in the hope that it will be useful,
17  but WITHOUT ANY WARRANTY; without even the implied warranty of
18  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19  Lesser General Public License for more details.
20 
21  You should have received a copy of the GNU Lesser General Public License
22  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
23 
24 *******************************************************************************
25 */
26 //@HEADER
27 
28 /*!
29  @file
30  @brief This file contains a class implementing a preconditioner for the Fluid-Structure Interaction system
31  @author Paolo Crosetto <crosetto@iacspc70.epfl.ch>
32  @date 08 Jun 2010
33 
34  We call the monolithic GCE matrix (with the fluid snd structure blocks C and N, couplings B and D):
35  \f$
36  A=\left(\begin{array}{cc}
37  C&B\\
38  D&N
39  \end{array}\right)\f$.
40 
41  The preconditioner is obtained from a block Gauss-Seidel preconditioner
42  \f$P_{DN}=
43  \left(\begin{array}{cc}
44  C&B\\
45  0&N
46  \end{array}\right)\f$,
47  decomposing it into two factors
48  \f$P_{DN}=P_1P_2
49  \left(
50  \begin{array}{cc}
51  I&0\\
52  0&N
53  \end{array}\right)
54  \left(\begin{array}{cc}
55  C&B\\
56  0&I
57  \end{array}\right)\f$ and applying a preconditioning strategy (algebraic additive Schwarz \f$P_{AS}\f$)
58  to each factor, so that \f$ P^{-1}=(P_{AS}(P_2))^{-1}(P_{AS}(P_1))^{-1}\f$.
59 
60  NOTE: this class is used also in the geometry implicit case, with an additional factor for the mesh motion, which is
61  coupled to the solid block using the default coupling method MonolithicBlockComposed::coupler. In that case the
62  preconditioner is decomposed in three factors, the fluid and structure ones being the same as for the GCE case.
63  NOTE2: this class is also the base class for other types of preconditioners, like ComposedDN2, ComposedDND. In fact for
64  instance it is used as F-S block in the preconditioners for the GI matrix in FSIFSIMonolithicGI
65  */
66 
67 #ifndef COMPOSEDDN_H
68 #define COMPOSEDDN_H 1
69 
70 #include <lifev/core/LifeV.hpp>
71 
72 #include <lifev/core/algorithm/PreconditionerComposed.hpp>
73 
74 #include <lifev/fsi/solver/MonolithicBlockComposed.hpp>
75 
76 namespace LifeV
77 {
78 
79 //! MonolithicBlockComposedDN - Short description of the class
80 /*!
81  @author Paolo Crosetto
82  @see \cite CrosettoEtAl2009
83 
84 
85  */
86 class MonolithicBlockComposedDN : public MonolithicBlockComposed
87 {
88 public:
89  typedef MonolithicBlockComposed super_Type;
90 
91  MonolithicBlockComposedDN ( const std::vector<Int>& flag, const std::vector<Int>& order) :
92  super_Type ( flag, order ),
93  M_blockPrecs()
94  {
95  }
96 
97  //! @name public methods
98  //@{
99 
100  //! sets the parameters related to M_blockPrecs from the data file
101  /*!
102  \param dataFile: GetPot data file
103  \param section: string identifying the section in the data file
104  */
105  void setDataFromGetPot ( const GetPot& dataFile,
106  const std::string& section );
107 
108  //! Solves the preconditioned linear system
109  /*!
110  Provided the linear solver and the right hand side this method computes the solution and returns it into
111  the result vector.
112  @param rhs right hand side of the linear system
113  @param result output result
114  @param linearSolver the linear system
115  */
116  int solveSystem ( const vector_Type& rhs, vector_Type& step, solverPtr_Type& linearSolver);
117 
118  //! Computes the coupling
119  /*!
120  computes all the coupling blocks specific for the chosen preconditioner. The coupling is handled
121  through an augmented formulation, introducing new variables (multipliers). Needs as input: the global map of the problem,
122  the FESpaces of the subproblems to be coupled with their offsets, a std::map holding the two numerations of the
123  interface between the two subproblems (the numeration can be different but the nodes must be matching in each
124  subdomain), an EpetraVector defined on the multipliers map containing the corresponding dof at the interface (NB: the multipliers map should be constructed from the second numeration in the std::map).
125 
126  In this case the coupling matrices are two:
127  \f$
128  C_1=
129  \left(
130  \begin{array}{cc}
131  I&0\\
132  0&0
133  \end{array}
134  \right)\f$
135  and
136  \f$
137  C_2=
138  \left(
139  \begin{array}{cc}
140  0&C\\
141  0&I
142  \end{array}
143  \right)\f$
144 
145  Note that the FESpaces and the offsets have to be set before calling this method.
146  @param map the map of the global problem
147  @param locDofMap std::map with the correspondence between the interface dofs for the two different maps in
148  the subproblems
149  @param numerationInterface vector containing the correspondence of the Lagrange multipliers with the interface dofs
150  */
151  virtual void coupler (mapPtr_Type& map,
152  const std::map<ID, ID>& locDofMap,
153  const vectorPtr_Type& numerationInterface,
154  const Real& timeStep,
155  const Real& coefficient,
156  const Real& rescaleFactor);
157 
158  //!pushes back the preconditioner for a block
159  /*!
160  In this case M_blockPrecs is of type PreconditionerComposed, thus this method calls PreconditionerComposed::push_back(...)
161  which computes the AAS preconditioner for the input matrix
162  \param Mat: input matrix
163  */
164  virtual void push_back_precs ( matrixPtr_Type& Mat);
165 
166  //! returns the true if the preconditioner has at leas one factor computed
167  bool set()
168  {
169  return M_blockPrecs->preconditionerCreated();
170  }
171 
172  /*! copies the shared_ptr to the communicator in the member M_comm and builds the empty ComposedPreconditioneronditioner
173  M_blockPrecs
174  */
175  void setComm ( std::shared_ptr<Epetra_Comm> comm )
176  {
177  M_comm = comm;
178  M_blockPrecs.reset ( new PreconditionerComposed (M_comm) );
179  }
180 
181  const std::vector<std::shared_ptr<Preconditioner> >& blockPrecs() const
182  {
183  return M_blockPrecs->composedPreconditionerPtr()->Operator();
184  }
185 
186  //@}
187  //!@name Factory Methods
188  //@{
189 
190  static MonolithicBlock* createComposedDN()
191  {
192  const Int order[] = { MonolithicBlockComposed::solid, MonolithicBlockComposed::fluid};
193  const Int couplingsDN[] = { 0, 7};
194  const std::vector<Int> couplingVectorDN (couplingsDN, couplingsDN + 2);
195  const std::vector<Int> orderVector (order, order + 2);
196  return new MonolithicBlockComposedDN (couplingVectorDN, orderVector);
197  }
198 
199 
200  static MonolithicBlock* createComposedDN2()
201  {
202  const Int order[] = { MonolithicBlockComposed::fluid, MonolithicBlockComposed::solid};
203  const Int couplingsDN2[] = { 8, 6};
204  const std::vector<Int> couplingVectorDN2 (couplingsDN2, couplingsDN2 + 2);
205  const std::vector<Int> orderVector (order, order + 2);
206  return new MonolithicBlockComposedDN (couplingVectorDN2, orderVector);
207  }
208 
209  static MonolithicBlock* createComposedDNGI()
210  {
211  //! Factorization in three:
212  /*!
213  - Solid: Neumann
214  | I | 0 | 0 | 0 |
215  |---+---+---+---|
216  | 0 | S | 0 | 0 |
217  |---+---+---+---|
218  | 0 | 0 | I | 0 |
219  |---+---+---+---|
220  | 0 | 0 | 0 | I |
221  - ALE: Dirichlet
222  | I | 0 | 0 | 0 |
223  |---+----+---+---|
224  | 0 | I | 0 | 0 |
225  |---+----+---+---|
226  | 0 | 0 | I | 0 |
227  |---+----+---+---|
228  | 0 | -C | 0 | H |
229  - Fluid: Dirichlet
230  | F | 0 | C | SD |
231  |---+-------+---+----|
232  | 0 | S | 0 | 0 |
233  |---+-------+---+----|
234  | C | -C/dt | 0 | 0 |
235  |---+-------+---+----|
236  | 0 | 0 | 0 | I |
237  */
238  const Int order[] = { MonolithicBlockComposed::solid, MonolithicBlockComposed::mesh, MonolithicBlockComposed::fluid };
239  const Int couplingsDNGI[] = { 0/*solid*/, 7/*fluid*/, 16/*ALE*/ };
240  const std::vector<Int> couplingVectorDNGI (couplingsDNGI, couplingsDNGI + 3);
241  const std::vector<Int> orderVector (order, order + 3);
242  return new MonolithicBlockComposedDN ( couplingVectorDNGI, orderVector );
243  }
244 
245 
246  static MonolithicBlock* createComposedDN2GI()
247  {
248  //! Factorization in three:
249  /*!
250  - Fluid: Dirichlet
251  | F | 0 | 0 | SD |
252  |---+---+---+----|
253  | 0 | I | 0 | 0 |
254  |---+---+---+----|
255  | C | 0 | I | 0 |
256  |---+---+---+----|
257  | 0 | 0 | 0 | I |
258  - Solid: Neumann
259  | I | 0 | 0 | 0 |
260  |---+-------+---+---|
261  | 0 | S | C | 0 |
262  |---+-------+---+---|
263  | 0 | -C/dt | 0 | 0 |
264  |---+-------+---+---|
265  | 0 | 0 | 0 | I |
266  - ALE: Dirichlet
267  | I | 0 | 0 | 0 |
268  |---+----+---+---|
269  | 0 | I | 0 | 0 |
270  |---+----+---+---|
271  | 0 | 0 | I | 0 |
272  |---+----+---+---|
273  | 0 | -C | 0 | H |
274  */
275  const Int order[] = { MonolithicBlockComposed::fluid, MonolithicBlockComposed::solid, MonolithicBlockComposed::mesh };
276  const Int couplingsDN2GI[] = { 8, 6, 16 };
277  const std::vector<Int> couplingVectorDN2GI (couplingsDN2GI, couplingsDN2GI + 3);
278  const std::vector<Int> orderVector (order, order + 3);
279  return new MonolithicBlockComposedDN ( couplingVectorDN2GI, orderVector );
280  }
281 
282  //@}
283 
284 protected:
285 
286  //! @name Protected Methods
287  //@{
288 
289  /*!
290  Replaces the preconditioner in M_blockPrecs with another one that is already constructed
291  */
292  virtual void replace_precs ( matrixPtr_Type& Mat, UInt position);
293  //@}
294 
295  //! @name Protected Members
296  //@{
297 
298  /*!
299  Pointer to an PreconditionerComposed object containing the preconditioners for each block
300  */
301  std::shared_ptr<PreconditionerComposed> M_blockPrecs;
302  //@}
303 
304 private:
305  // static bool reg;
306 
307 };
308 
309 } // Namespace LifeV
310 
311 #endif /* COMPOSEDDN_H */
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
double Real
Generic real data.
Definition: LifeV.hpp:175