LifeV
MonolithicBlockComposedDNND.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 File containing a calss for Neumann-Neumann preconditioners
31 
32  @author Paolo Crosetto <crosetto@iacspc70.epfl.ch>
33  @date 02 Jul 2010
34 
35  File containing a calss for composed preconditioner of the following type:
36  given the matrix \f$A=A_1A_2+A_3A_4f$, \f$A^{-1}\approxA_1^{-1}A_2^{-1}+A_3^{-1}A_4^{-1}\approx P_1^{-1}P_2^{-1}+P_3^{-1}P_4^{-1}\f$ then we compute the preconditioner
37  \f$P^{-1}=_4^{-1}P_3^{-1}+P_2^{-1}P_1^{-1}\f$.
38  */
39 
40 #ifndef COMPOSEDDNND_H
41 #define COMPOSEDDNND_H 1
42 
43 #include <lifev/core/LifeV.hpp>
44 #include <lifev/fsi/solver/MonolithicBlockComposedNN.hpp>
45 
46 namespace LifeV
47 {
48 
49 //! ComposedPrecDNND - Short description of the class
50 /*!
51  @author Paolo Crosetto
52 
53  Class implementing a composed preconditioner of the following type:
54  given the matrix \f$A= A_1 A_2+ A_3 A_4\approx P_1 P_2+ P_3 P_4\f$ then we compute the preconditioner
55  \f$P^{-1}=P_4^{-1}P_3^{-1}+P_2^{-1}P_1^{-1}\f$.
56  In particular in this case we use for \f$A_1\f$ and \f$A_4\f$ Dirichlet problems, for \f$A_2\f$ and \f$A_3\f$ Neumann problems. The form of the decomposition is the following:
57 \f$
58 A=A_{(1)}+A_{(2)}
59 =
60 \left(
61 \begin{array}{ccc}
62  2C& 0 & \mathcal I\\
63  0&2N& \mathcal I\\
64 0&0&I
65 \end{array}
66 \right)
67 +
68 \left(
69 \begin{array}{ccc}
70 2C& 0 & 0\\
71 0 &2N&0\\
72 \tilde\Delta_1&\tilde\Delta_2&-I
73 \end{array}
74 \right).
75 \f$
76 So a block factorization for these two matrices can easily be computed,
77 in fact
78 \f$
79 A_{(1)}
80 =
81 \left(
82 \begin{array}{ccc}
83  2C& 0 & \mathcal I\\
84 0 &I&0\\
85 2\tilde\Delta_1&0&0
86 \end{array}
87 \right)
88 \left(
89 \begin{array}{ccc}
90 I& 0 & 0 \\
91 0 &2 N&2\mathcal I\\
92 0&0&I
93 \end{array}
94 \right)
95 \f$
96 while
97 \f$
98 A_{(2)}
99 =
100 \left(
101 \begin{array}{ccc}
102  2C& 0 & 0\\
103 0 &I&0\\
104 2\tilde\Delta_1&0&I
105 \end{array}
106 \right)
107 \left(
108 \begin{array}{ccc}
109 I& 0 & 0 \\
110 0 &2 N&2\mathcal I\\
111 0&\tilde\Delta_2&0
112 \end{array}
113 \right).
114 \f$
115  */
116 class MonolithicBlockComposedDNND : public MonolithicBlockComposedNN
117 {
118 public:
119 
120  typedef MonolithicBlockComposedNN super_Type;
121 
122  //! @name Constructor and Destructor
123  //@{
124 
125  MonolithicBlockComposedDNND ( const std::vector<Int>& flag, const std::vector<Int>& order ) :
126  super_Type ( flag, order )
127  {
128  }
129 
130  ~MonolithicBlockComposedDNND()
131  {}
132 
133  //@}
134 
135  //! @name Public Methods
136  //@{
137 
138  //! Computes the coupling
139  /*!
140  computes all the coupling blocks specific for the chosen preconditioner. The coupling is handled
141  through an augmented formulation, introducing new variables (multipliers). Needs as input: the global map of the problem,
142  the FESpaces of the subproblems to be coupled with their offsets, a std::map holding the two numerations of the
143  interface between the two subproblems (the numeration can be different but the nodes must be matching in each
144  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).
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  @param timestep the timestep chosen
151  */
152  void coupler (mapPtr_Type& map,
153  const std::map<ID, ID>& locDofMap,
154  const vectorPtr_Type& numerationInterface,
155  const Real& timeStep,
156  const Real& coefficient,
157  const Real& rescaleFactor);
158 
159  //@}
160  //!@name Factory Method
161  //@{
162 
163  static MonolithicBlock* createComposedDNND()
164  {
165  const Int couplingsDNND[] = { 8, 4, 2, 8, 1, 2 };
166  const Int order[] = { MonolithicBlockComposedNN::fluid1, MonolithicBlockComposedNN::solid1, MonolithicBlockComposedNN::fluid2, MonolithicBlockComposedNN::solid2};
167  const std::vector<Int> couplingVectorDNND (couplingsDNND, couplingsDNND + 6);
168  const std::vector<Int> orderVector (order, order + 4);
169  return new MonolithicBlockComposedDNND (couplingVectorDNND, orderVector);
170  }
171 
172  //@}
173 
174 private:
175 
176  //static bool reg;
177 };
178 
179 } // Namespace LifeV
180 
181 #endif /* COMPOSEDDNND_H */
void assignFunction(bcBase_Type &base)
Assign the function to the base of the BCHandler.
int32_type Int
Generic integer data.
Definition: LifeV.hpp:188
double Real
Generic real data.
Definition: LifeV.hpp:175