DarcySolverNonLinear.hpp 26.7 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
//@HEADER
/*
*******************************************************************************

   Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
   Copyright (C) 2010 EPFL, Politecnico di Milano, Emory UNiversity

   This file is part of the LifeV library

   LifeV is free software; you can redistribute it and/or modify
   it under the terms of the GNU Lesser General Public License as published by
   the Free Software Foundation; either version 3 of the License, or
   (at your option) any later version.

   LifeV is distributed in the hope that it will be useful,
   but WITHOUT ANY WARRANTY; without even the implied warranty of
   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
   Lesser General Public License for more details.

   You should have received a copy of the GNU Lesser General Public
   License along with this library; if not, see <http://www.gnu.org/licenses/>
22
23
24
25
26
27
28
29
30

 This file is part of the LifeV library

 Author(s): A. Fumagalli  <alessio.fumagalli@mail.polimi.it>
      Date: 2010-05-09

 Copyright (C) 2001-2006 EPFL, Politecnico di Milano, INRIA
 Copyright (C) 2006-2010 Politecnico di Milano

31
*******************************************************************************
32
33

*/
34
//@HEADER
35

36
37
/*!
  @file
38
  @brief This file contains a non-linear permeability term Darcy equation solver class
39
40
41
42
43
44

  @date 05/2010
  @author A. Fumagalli <alessio.fumagalli@mail.polimi.it>

  @contributor M. Kern <michel.kern@inria.fr>
  @maintainer M. Kern <michel.kern@inria.fr>
45
46
*/

47
48
#ifndef _DARCYSOLVERNONLINEAR_H_
#define _DARCYSOLVERNONLINEAR_H_ 1
49

50
#include <lifev/darcy/solver/DarcySolver.hpp>
51
52
53
54
55
56
57
58
59
60
61

// Local namespace to store the default function for start up the fixed point scheme and the default non-linear permeability tensor.
namespace
{

/*! @class Default Darcy start up function for the fixed point scheme. Needed in DarcySolverNonLinear class.
  This class implement the default start up function for the fixed point scheme in the the non-linear Darcy problem,
  it is the zero function.
*/
struct DarcyDefaultStartUpFunction
{
62
    LifeV::Real operator()( const LifeV::Real&, const LifeV::Real&,
63
64
                const LifeV::Real&, const LifeV::Real&,
                const LifeV::UInt&) const
65
    {
Tiziano Passerini's avatar
Tiziano Passerini committed
66
        return static_cast<LifeV::Real>( 0 );
67
68
69
70
71
72
73
74
    }
};

}

// LifeV namespace.
namespace LifeV
{
75
76
//!  @class DarcySolverNonLinear This class implements a non-linear Darcy solver

77
/*!
78
  @author A. Fumagalli <alessio.fumagalli@mail.polimi.it>
79
  @see For applications related to two-phase flow see \cite Fumagalli2011a
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250

  This class implements a non-linear Darcy solver with fixed point scheme to handle the non-linearity.
  <br>
  The classical non-linear Darcy formulation is a couple of differential equations of first order with
  the unknowns \f$ p \in C^1 (\Omega ) \f$, being the pressure or the primal unknown,
  and \f$ \sigma \in (C^1( \Omega ) )^n \f$, being the Darcy velocity or the flux or the dual unknown,
  such that
  \f[
  \left\{
  \begin{array}{l l l}
  \Lambda^{-1}(p) \sigma + \nabla p = 0 & \mathrm{in} & \Omega\,,  \vspace{0.2cm} \\
  \nabla \cdot \sigma - f = 0        & \mathrm{in} & \Omega\,,  \vspace{0.2cm} \\
  p = g_D                            & \mathrm{on} & \Gamma_D\,,\vspace{0.2cm} \\
  \sigma \cdot n + h p = g_R         & \mathrm{on} & \Gamma_R\,.
  \end{array}
  \right.
  \f]
  Where \f$ \Lambda(p) \f$ is the permeability tensor that depends on \f$ p \f$, \f$ f \f$ is the source term,
  \f$ \Gamma_D \f$ is the subset  of the boundary of \f$ \Omega \f$ with Dirichlet boundary conditions with datum
  \f$ g_D \f$ and \f$ \Gamma_R \f$  is the part of the boundary of \f$ \Omega \f$ with Robin, or Neumann, boundary
  conditions with data \f$ h \f$ and \f$ g_R \f$. We suppose that \f$ \partial \Omega = \Gamma_D \cup \Gamma_R \f$ and
  \f$ \Gamma_D \cap \Gamma_R = \emptyset \f$.
  <br>
  Using the hybridization procedure, and introducing a new variable, we may split the problem from
  the entire domain to problems in each elements of the triangulation \f$ \mathcal{T}_h \f$, then we may write
  the weak formulation for the dual problem.
  The new variable is the hybrid unknown \f$ \lambda \f$, e.g. "the trace" of pressure of the primal unknown,
  it is the Lagrange multipliers forcing the continuity of the flux across each faces or edges in the triangulation.
  We introduce the following functional spaces, the first is the space of the primal variable, the third for
  the dual variable, the fourth for the hybrid variable.
  \f[
  \begin{array}{l}
  V = L^2 (\Omega ) \,,\vspace{0.2cm} \\
  H(div, K) = \left\{ \tau \in ( L^2(K))^n : \nabla \cdot \tau \in L^2(K)\right\}\,, \vspace{0.2cm}\\
  Z = \left\{ \tau \in L^2(\Omega) : \tau\vert_K \in H(div, K) \, \forall K \in  \mathcal{T}_h \right\}\,, \vspace{0.2cm}\\
  \Lambda = \left\{ \lambda \in \prod_{K \in \mathcal{T}_h} H^{1/2} (\partial K): \lambda_K = \lambda_{K'} \,\, \mathrm{on} \,\, e_{K-K'} \, \forall K \in \mathcal{T}_h,\, \lambda = g_D \,\, \mathrm{on} \,\, \Gamma_D \right\}\,.
  \end{array}
  \f]
  Introducing the following bilinear forms, functionals and operator
  \f[
  \begin{array}{l l}
  a(\sigma, \tau, p) = \sum_{K \in \mathcal{T}_h}\int_K \Lambda^{-1}(p) \sigma \cdot \tau \,,  &
  b(p, \tau) = -\sum_{K \in \mathcal{T}_h} \int_K p \nabla \cdot \tau\,, \vspace{0.2cm}\\
  c(\lambda, \tau) = \sum_{K \in \mathcal{T}_h} \int_{\partial K} \lambda \tau \cdot n\,,&
  h(\lambda, \mu) = \sum_{K \in \mathcal{T}_h} \int_{\partial K \cap \Gamma_R} h \mu \lambda \,,\vspace{0.2cm}\\
  F(v) = \sum_{K \in \mathcal{T}_h} \int_K f v\,,&
  G(\mu) =\sum_{K \in \mathcal{T}_h} \int_{\partial K \cap \Gamma_R} g \mu\,,
  \end{array}
  \f]
  we obtain the Darcy problem in the weak form: find \f$ (\sigma, \, p, \, \lambda) \in Z \times V \times \Lambda \f$ such that
  \f[
  \left\{
  \begin{array}{l l}
  a(\sigma, \tau, p) + b(p, \tau) + c(\lambda, \tau) = 0\,,  & \forall \tau \in Z \,,\vspace{0.2cm}\\
  b(v, \sigma) = - F(v)\,,                                  & \forall v \in V \,,\vspace{0.2cm}\\
  c(\mu, \sigma) + h(\lambda, \mu) = G(\mu) \,,           & \forall \mu \in \Lambda\,.
  \end{array}
  \right.
  \f]
  At discrete level we introduce the polynomial space, of degree \f$ r \f$, that approximate the finite dimensional
  spaces introduced above \f$ V_h \subset V \f$, \f$ Z_h \subset Z \f$ and \f$\Lambda_h \subset \Lambda \f$
  \f[
  \begin{array}{l}
  V_h = \left\{ v_h \in V: v_h|_K \in P_r (K)\, \forall K \in \mathcal{T}_h \right\}\,, \vspace{0.2cm}\\
  Z_h = \left\{ \tau_h \in Z: \tau_h|_K \in RT_r(K) \, \forall K \in \mathcal{T}_h \right\} \,, \vspace{0.2cm} \\
  \Lambda_h = \left\{ \lambda_h \in \Lambda: \lambda_h|_{\partial K} \in R_r( \partial K ) \, \forall K \in \mathcal{T}_h\right\}\,,
  \end{array}
  \f]
  where \f$ P_r(K) \f$ is the space of polynomial of degree \f$ r \f$ in the element \f$ K \f$, \f$ RT_r(K) \f$ is the space of
  polynomial of Raviart-Thomas of degrees \f$ r \f$ in the element \f$ K \f$ and \f$ R_r(\partial K) \f$ is the space of polynomial of
  degree \f$ r \f$ definite on each of the boundary of the element \f$ K \f$ and discontinuous from one edge to the other.
  <BR>
  The finite dimensional problem is: find \f$ (\sigma_h,\,, p_h, \, \lambda_h) \in Z_h \times V_h \times \Lambda_h \f$ such that
  \f[
  \left\{
  \begin{array}{l l}
  a(\sigma_h, \tau_h, p_h) + b(p_h, \tau_h) + c(\lambda_h, \tau_h) = 0\,,  & \forall \tau_h \in Z_h \,,\vspace{0.2cm}\\
  b(v_h, \sigma_h) = - F(v_h)\,,                                  & \forall v_h \in V_h \,,\vspace{0.2cm}\\
  c(\mu_h, \sigma_h) + h(\lambda_h, \mu_h) = G(\mu_h) \,,           & \forall \mu_h \in \Lambda_h\,.
  \end{array}
  \right.
  \f]
  To solve the non-linearity we use a fixed point scheme based on the relative difference between two consecutive iterations
  of the primal variable. We start from the, user defined, function \f$ p^0 \f$ and solve the linearized problem for \f$ n \geq 1 \f$:
  find \f$ (\sigma_h^n,\, p_h^n, \, \lambda_h^n) \in Z_h \times V_h \times \Lambda_h \f$ such that
  \f[
  \left\{
  \begin{array}{l l}
  a(\sigma_h^n, \tau_h, p_h^{n-1}) + b(p_h^n, \tau_h) + c(\lambda_h^n, \tau_h) = 0\,,  & \forall \tau_h \in Z_h \,,\vspace{0.2cm}\\
  b(v_h, \sigma_h^n) = - F(v_h)\,,                                  & \forall v_h \in V_h \,,\vspace{0.2cm}\\
  c(\mu_h, \sigma_h^n) + h(\lambda_h^n, \mu_h) = G(\mu_h) \,,           & \forall \mu_h \in \Lambda_h\,.
  \end{array}
  \right.
  \f]
  At each iteration, to solve the problem, we use the static condensation procedure, i.e. the unknowns in the discrete
  weak system are not independent and \f$ p_K^n \f$, \f$\sigma_K^n \f$ may be written in function of
  \f$ \lambda_K^n \f$ alone. We introduce the following local matrices
  \f[
  \begin{array}{l l l}
  \left[ A \left( p_K^{n-1}\right) \right]_{ij} =   \int_K \Lambda^{-1}\left( p^{n-1} \right) \psi_j \cdot \psi_i \,, &
  \left[ B \right]_{ij} = - \int_K \phi_j \nabla \cdot \psi_i \,, &
  \left[ C \right]_{ij} =   \int_{\partial K} \xi_i \psi_j \cdot n \,, \vspace{0.2cm} \\
  \left[ H \right]_{ij} =   \int_{\partial K \cap \Gamma_R} h \xi_i \xi_j\,, &
  \left[ F \right]_{j}  =   \int_K f \phi_j\,, &
  \left[ G \right]_{j}  =   \int_{\partial K \cap \Gamma_R } g \xi_j\,,
  \end{array}
  \f]
  where we avoid to write the dependence on the triangle \f$ K \f$ in all the matrices and vectors.
  <br>
  The local matrix formulation of the finite dimensional problem, at each iteration, is
  \f[
  \left\{
  \begin{array}{l}
  A \left( p_K^{n-1} \right) \sigma_K^n + B p_K^n + C \lambda_K^n = 0\,, \vspace{0.2cm} \\
  B^T \sigma_K^n = -F \,,                    \vspace{0.2cm}\\
  C^T \sigma_K^n + H \lambda_K^n = G\,.
  \end{array}
  \right.
  \f]
  Or alternatively
  \f[
  \begin{array}{l l l}
  \left[
  \begin{array}{c c c}
  A \left( p_K^{n-1} \right) & B & C \vspace{0.2cm} \\
  B^T                        & 0 & 0 \vspace{0.2cm} \\
  C^T                        & 0 & H
  \end{array}
  \right] \, \cdot &
  \left[
  \begin{array}{c}
  \sigma_K^n  \vspace{0.2cm}\\
  p_K^n       \vspace{0.2cm}\\
  \lambda_K^n
  \end{array}
  \right] \, &
  =
  \left[
  \begin{array}{c}
  0 \vspace{0.2cm}\\
  -F \vspace{0.2cm}\\
  G
  \end{array}
  \right]\,.
  \end{array}
  \f]
  Introducing the local hybrid matrix and local hybrid right hand side
  \f[
  \begin{array}{l}
  L_K =  -C^T A^{-1}\left( p_k^{n-1}\right) C + C^T A^{-1} \left( p_K^{n-1} \right) B \left[ B^T A^{-1}\left( p_K^{n-1}\right) B \right]^{-1} B^T A^{-1}\left( p_K^{n-1} \right) C + H \,, \vspace{0.2cm} \\
  r_K = G + C^T A^{-1}\left( p_K^{n-1} \right) B \left[ B^T A^{-1}\left( p_K^{n-1} \right) B \right]^{-1} F\,,
  \end{array}
  \f]
  Imposing that at each edge or face the hybrid unknown is single value we obtain a linear system for the hybrid unknown
  \f[
  L^n \lambda^n = r \,.
  \f]
  We recover the primal and dual variable as a post-process from the hybrid variable at element level, so we have
  \f[
  \begin{array}{l}
  p_K^n = \left[ B^T A^{-1}\left( p_K^{n-1} \right) B \right]^{-1} \left[ F -  B^T A^{-1}\left( p_K^{n-1}\right) C \lambda_K^n \right]\,, \vspace{0.2cm} \\
  \sigma_K^n = -A^{-1}\left(p_K^{n-1} \right) \left( B p_K^n + C \lambda_K^n \right) \,.
  \end{array}
  \f]
  @note In the code we do not use the matrix \f$ H \f$ and the vector \f$ G \f$, because all the boundary
  conditions are imposed via BCHandler class.
  @todo Insert any scientific publications that use this solver.
  @todo Attention! We have the hypothesis that we use P0 elements for the primal unknown. Change this in a future!
  @todo Add criteria to ensure convergence of the fixed point method.
*/
template< typename Mesh,
251
typename SolverType = LifeV::SolverAztecOO >
252
class DarcySolverNonLinear :
253
254
255
256
257
        virtual public DarcySolver<Mesh, SolverType>
{

public:

258
    //! @name Public Types
259
260
261
262
    //@{

    typedef DarcySolver<Mesh, SolverType> DarcySolverPolicies;

fumagalli's avatar
fumagalli committed
263
    typedef typename DarcySolverPolicies::function_Type      function_Type;
264

265
    typedef typename DarcySolverPolicies::permeability_Type  permeability_Type;
266

267
    typedef typename DarcySolverPolicies::data_Type          data_Type;
268

269
270
    typedef typename DarcySolverPolicies::bchandler_raw_Type bchandler_raw_Type;
    typedef typename DarcySolverPolicies::bchandler_Type     bchandler_Type;
271

272
273
    typedef typename DarcySolverPolicies::matrix_Type        matrix_Type;
    typedef typename DarcySolverPolicies::matrixPtr_Type     matrixPtr_Type;
274

275
276
    typedef typename DarcySolverPolicies::vector_Type        vector_Type;
    typedef typename DarcySolverPolicies::vectorPtr_Type     vectorPtr_Type;
277

278
279
    typedef typename DarcySolverPolicies::comm_Type          comm_Type;
    typedef typename DarcySolverPolicies::commPtr_Type       commPtr_Type;
280
281
282
283
284
285
286
287
288
289
290
291
292
293

    //@}

    //! @name Constructors & destructor
    //@{

    /*!
      Full constructor for the class.
      @param dataFile Data for the problem.
      @param primal_FESpace Primal finite element space.
      @param dual_FESpace Dual element space.
      @param hybrid_FESpace Hybrid finite element space.
      @param VdotN_FESpace Dual basis function dot outward unit normal at each face (3D) or edge (2D) finite element space.
      @param bcHandler Boundary conditions for the problem.
294
      @param comm Shared pointer of the Epetra communicator.
295
    */
296
    DarcySolverNonLinear ( const data_Type&          dataFile,
297
298
299
300
                           FESpace<Mesh, MapEpetra>& primal_FESpace,
                           FESpace<Mesh, MapEpetra>& dual_FESpace,
                           FESpace<Mesh, MapEpetra>& hybrid_FESpace,
                           FESpace<Mesh, MapEpetra>& VdotN_FESpace,
301
                           bchandler_raw_Type&       bcHandler,
Michel Kern's avatar
Michel Kern committed
302
                           const commPtr_Type&             comm );
303
304
305
306
307
308
309
310
311

    /*!
      Constructor for the class without the definition of the boundary handler.
      @param dataFile Data for the problem.
      @param primal_FESpace Primal finite element space.
      @param dual_FESpace Dual finite element space.
      @param hybrid_FESpace Hybrid finite element space.
      @param VdotN_FESpace Dual basis function dot outward unit normal at each face (3D) or edge (2D) finite element space.
      @param bcHandler Boundary conditions for the problem.
312
      @param comm Shared pointer of the Epetra communicator.
313
    */
314
    DarcySolverNonLinear ( const data_Type&          dataFile,
315
316
317
318
                           FESpace<Mesh, MapEpetra>& primal_FESpace,
                           FESpace<Mesh, MapEpetra>& dual_FESpace,
                           FESpace<Mesh, MapEpetra>& hybrid_FESpace,
                           FESpace<Mesh, MapEpetra>& VdotN_FESpace,
Michel Kern's avatar
Michel Kern committed
319
                           const commPtr_Type&             comm );
320
321

    //! Virtual destructor.
322
    virtual ~DarcySolverNonLinear ();
323
324
325

    //@}

326
    //! @name Methods
327
328
    //@{

329
    //!  Set up the linear solver and the preconditioner for the linear system.
330
    virtual void setup ();
331

332
333
334
335
336
337
338
339
340
    //!  Perform the fixed point scheme.
    void fixedPointScheme ();

    //@}

    //! @name Set methos
    //@{

    //! Initial guess for fixed point iteration
341
342
343
344
    /*!
      Set the function for the first iteration for the fixed point method. The default is the zero function.
      @param primalZeroIteration The function for the first iteration.
    */
fumagalli's avatar
fumagalli committed
345
    void setPrimalZeroIteration ( const function_Type& primalZeroIteration );
fumagalli's avatar
fumagalli committed
346

347
    //! Set the inverse of diffusion tensor,
348
    /*!
349
      The default inverse of permeability is the identity matrix.
350
351
      @param invPerm Inverse of the permeability tensor for the problem.
    */
Michel Kern's avatar
Michel Kern committed
352
    void setInversePermeability ( const permeability_Type& invPerm )
353
354
355
356
357
358
    {

        // Call the standard set inverse permeability
        DarcySolver<Mesh, SolverType>::setInversePermeability( invPerm );

        // Set the dependence of the previous solution in the permeability
Michel Kern's avatar
Michel Kern committed
359
        this->M_inversePermeability->setField( primalPreviousIteration() );
360
361
    }

362
363
364
365
366
367
368
369
370
    //@}

    //! @name Get methods
    //@{

    /*!
      Returns the number of iteration of the fixed point scheme.
      @return Final number of iterations of the fixed point method as a constant UInt.
    */
Michel Kern's avatar
Michel Kern committed
371
    UInt fixedPointNumIteration ()
372
    {
fumagalli's avatar
fumagalli committed
373
374
375
        return M_fixedPointNumIteration;
    }

376
    //!  Returns the residual between two iterations of the fixed point scheme.
377
378
379
    /*!
      @return Final residual of the fixed point method as a constant Real.
    */
Michel Kern's avatar
Michel Kern committed
380
    Real fixedPointResidual ()
fumagalli's avatar
fumagalli committed
381
382
383
384
    {
        return M_fixedPointResidual;
    }

Michel Kern's avatar
Michel Kern committed
385
386
387
388
    //! Returns fixed point tolerance
    /*!
      @return M_fixedPointTolerance
     */
389
    const Real& fixedPointTolerance () const
Michel Kern's avatar
Michel Kern committed
390
    {
fumagalli's avatar
fumagalli committed
391
        return M_fixedPointTolerance;
Michel Kern's avatar
Michel Kern committed
392
    }
fumagalli's avatar
fumagalli committed
393
394
395
396
397
398

    //! Returns fixed point tolerance
    /*!
      @return M_fixedPointTolerance
    */
    Real& fixedPointTolerance ()
Michel Kern's avatar
Michel Kern committed
399
    {
fumagalli's avatar
fumagalli committed
400
        return M_fixedPointTolerance;
Michel Kern's avatar
Michel Kern committed
401
    }
402

Michel Kern's avatar
Michel Kern committed
403
404
405
406
    //! Returns maximum number of fixed point iterations allowed
    /*!
      @return max possible number of fixed point iterations
    */
407
    const UInt& fixedPointMaxIteration () const
Michel Kern's avatar
Michel Kern committed
408
    {
fumagalli's avatar
fumagalli committed
409
        return M_fixedPointMaxIteration;
410
    }
Michel Kern's avatar
Michel Kern committed
411

fumagalli's avatar
fumagalli committed
412
413
414
415
416
    //! Returns maximum number of fixed point iterations allowed
    /*!
      @return max possible number of fixed point iterations
    */
    UInt& fixedPointMaxIteration ()
Michel Kern's avatar
Michel Kern committed
417
    {
fumagalli's avatar
fumagalli committed
418
        return M_fixedPointMaxIteration;
419
    }
Michel Kern's avatar
Michel Kern committed
420

421
    //!  Returns the pointer of the primal solution vector at previous step.
422
    /*!
423
      @return Constant vector_Type reference of the primal solution at previous step.
424
    */
Michel Kern's avatar
Michel Kern committed
425
426
427
428
    const vectorPtr_Type& primalPreviousIteration () const
    {
        return M_primalPreviousIteration;
    }
429

Michel Kern's avatar
Michel Kern committed
430
          vectorPtr_Type& primalPreviousIteration ()
431
432
433
434
    {
        return M_primalPreviousIteration;
    }

435
436
437
438
    //@}

protected:

439
440
441
442
443
    // Methods
    //! @name Protected Methods
    //@{

    //! Update all problem variables
444
445
446
447
448
    /*!
      Update all the variables of the problem before the construction of
      the global hybrid matrix, e.g. reset the global hybrid matrix.
      It is principally used for a time dependent derived class.
    */
449
    virtual void updateVariables ();
450
451
452

    //@}

Michel Kern's avatar
Michel Kern committed
453
454
455
456
457
458
459
460
461
    //! @name Set methods
    //@{

    /*!
      Set fixed point tolerance
      @param tol requested tolerance
    */
    void setFixedPointTolerance ( const Real& tol)
    {
fumagalli's avatar
fumagalli committed
462
        M_fixedPointTolerance = tol;
Michel Kern's avatar
Michel Kern committed
463
464
465
466
467
468
469
470
    }

    /*!
      Set maximum number of fixed point iterations permitted
      @param maxit requested maximum
     */
    void setFixedPointMaxIteration ( const UInt& maxit)
    {
fumagalli's avatar
fumagalli committed
471
        M_fixedPointMaxIteration = maxit;
Michel Kern's avatar
Michel Kern committed
472
473
474
475
476
477
478
479
480
    }

    //@}

    //! @name protected variables
    //@{

    //! Primal solution at previous iteration step.
    vectorPtr_Type M_primalPreviousIteration;
481

Michel Kern's avatar
Michel Kern committed
482
483
484
    //@}

private:
485

486
487
488
489
490
    // Non-linear stuff.
    //! @name Non-linear stuff
    //@{

    //! The maximum number of iterations for the fixed point method.
fumagalli's avatar
fumagalli committed
491
    UInt           M_fixedPointMaxIteration;
492
493

    //! The current iterations for the fixed point method.
fumagalli's avatar
fumagalli committed
494
    UInt           M_fixedPointNumIteration;
495
496

    //! Tollerance for the stopping criteria for the fixed point method.
fumagalli's avatar
fumagalli committed
497
    Real           M_fixedPointTolerance;
498
499

    //! The residual between two iteration of the fixed point scheme.
fumagalli's avatar
fumagalli committed
500
    Real           M_fixedPointResidual;
501
502

    //! Primal solution at zero time step.
fumagalli's avatar
fumagalli committed
503
    function_Type  M_primalZeroIteration;
504
505
506
507
508
509
510
511
512

    //@}

}; // class DarcySolverNonLinear

//
// IMPLEMENTATION
//

513
514
515
516
// ===================================================
// Constructors & Destructor
// ===================================================

517
518
519
// Complete constructor.
template<typename Mesh, typename SolverType>
DarcySolverNonLinear<Mesh, SolverType>::
520
DarcySolverNonLinear ( const data_Type&           dataFile,
521
522
523
524
                       FESpace<Mesh, MapEpetra>&  primal_FESpace,
                       FESpace<Mesh, MapEpetra>&  dual_FESpace,
                       FESpace<Mesh, MapEpetra>&  hybrid_FESpace,
                       FESpace<Mesh, MapEpetra>&  VdotN_FESpace,
525
                       bchandler_raw_Type&        bcHandler,
Michel Kern's avatar
Michel Kern committed
526
                       const commPtr_Type&              comm ):
Tiziano Passerini's avatar
Tiziano Passerini committed
527
528
529
        // Standard Darcy solver constructor.
        DarcySolver<Mesh, SolverType>::DarcySolver( dataFile, primal_FESpace, dual_FESpace, hybrid_FESpace, VdotN_FESpace, bcHandler, comm),
        // Non-linear stuff.
fumagalli's avatar
fumagalli committed
530
        M_primalPreviousIteration       ( new vector_Type ( this->M_primal_FESpace.map() ) ),
fumagalli's avatar
fumagalli committed
531
532
533
534
        M_fixedPointMaxIteration        ( static_cast<UInt>(10) ),
        M_fixedPointNumIteration        ( static_cast<UInt>(0) ),
        M_fixedPointTolerance           ( static_cast<Real>(1.e-8) ),
        M_fixedPointResidual            ( M_fixedPointTolerance + static_cast<Real>(1) ),
Tiziano Passerini's avatar
Tiziano Passerini committed
535
        M_primalZeroIteration           ( DarcyDefaultStartUpFunction() )
536
537
538
539

{

    // Interpolate the primal initial value on the default initial value function, for the fixed point scheme.
540
541
542
    this->M_primal_FESpace.interpolate( M_primalZeroIteration,
                                        *(this->M_primal),
                                        static_cast<Real>(0) );
543
544
545
546
547
548

} // Constructor

// Constructor without boundary condition handler.
template<typename Mesh, typename SolverType>
DarcySolverNonLinear<Mesh, SolverType>::
549
DarcySolverNonLinear ( const data_Type&           dataFile,
550
551
552
553
                       FESpace<Mesh, MapEpetra>&  primal_FESpace,
                       FESpace<Mesh, MapEpetra>&  dual_FESpace,
                       FESpace<Mesh, MapEpetra>&  hybrid_FESpace,
                       FESpace<Mesh, MapEpetra>&  VdotN_FESpace,
Michel Kern's avatar
Michel Kern committed
554
                       const commPtr_Type&              comm ):
Tiziano Passerini's avatar
Tiziano Passerini committed
555
556
557
        // Standard Darcy solver constructor.
        DarcySolver<Mesh, SolverType>::DarcySolver( dataFile, primal_FESpace, dual_FESpace, hybrid_FESpace, VdotN_FESpace, comm),
        // Non-linear stuff.
fumagalli's avatar
fumagalli committed
558
        M_primalPreviousIteration       ( new vector_Type ( this->M_primal_FESpace.map() ) ),
fumagalli's avatar
fumagalli committed
559
560
561
562
        M_fixedPointMaxIteration        ( static_cast<UInt>(10) ),
        M_fixedPointNumIteration        ( static_cast<UInt>(0) ),
        M_fixedPointTolerance           ( static_cast<Real>(1.e-8) ),
        M_fixedPointResidual            ( M_fixedPointTolerance + static_cast<Real>(1) ),
Tiziano Passerini's avatar
Tiziano Passerini committed
563
        M_primalZeroIteration           ( DarcyDefaultStartUpFunction() )
564
565
566
{

    // Interpolate the primal initial value on the default initial value function, for the fixed point scheme.
567
568
569
    this->M_primal_FESpace.interpolate( M_primalZeroIteration,
                                        *(this->M_primal),
                                        static_cast<Real>(0) );
570
571
572
573
574
575

} // Constructor

// Virtual destructor.
template<typename Mesh, typename SolverType>
DarcySolverNonLinear<Mesh, SolverType>::
576
~DarcySolverNonLinear ( void )
577
578
579
580
{

} // Destructor

581
582
583
584
// ===================================================
// Public methods
// ===================================================

585
586
587
588
// Set up the linear solver and the preconditioner.
template<typename Mesh, typename SolverType>
void
DarcySolverNonLinear<Mesh, SolverType>::
589
setup ()
590
591
592
593
594
595
596
597
{

    GetPot dataFile( *(this->M_data.dataFile()) );

    // Call the DarcySolver setup method for setting up the linear solver.
    DarcySolver<Mesh, SolverType>::setup();

    // Set the maximum number of iteration for the fixed point iteration scheme.
Michel Kern's avatar
Michel Kern committed
598
    setFixedPointMaxIteration(static_cast<UInt>( dataFile( ( this->M_data.section()
599
                                 + "/non-linear/fixed_point_iteration" ).data(), 10 ) ));
600
601

    // Set the tollerance for the fixed point iteration scheme.
Michel Kern's avatar
Michel Kern committed
602
    setFixedPointTolerance(dataFile( ( this->M_data.section() + "/non-linear/fixed_point_toll" ).data(), 1.e-8 ));
603
604
605
606
607
608
609

} // setup

// Fixed point scheme.
template <typename Mesh, typename SolverType>
void
DarcySolverNonLinear<Mesh, SolverType>::
610
fixedPointScheme ()
611
612
613
{

    // Current iteration.
fumagalli's avatar
fumagalli committed
614
    M_fixedPointNumIteration = static_cast<UInt>(0);
615
616

    // Error between two iterations, it is the relative error between two step of the primal vector
Michel Kern's avatar
Michel Kern committed
617
    M_fixedPointResidual =  fixedPointTolerance() + 1;
618
619
620

    /* A loop for the fixed point scheme, with exit condition based on stagnate of the
       primal variable and the maximum iteration. */
621
    while (    fixedPointResidual() > fixedPointTolerance()
622
           && fixedPointNumIteration() < fixedPointMaxIteration() )
623
624
    {
        // Increment the iteration number.
fumagalli's avatar
fumagalli committed
625
        ++M_fixedPointNumIteration;
626
627
628
629
630
631
632
633
634
635
636

        // Build the linear system and the right hand side.
        this->buildSystem();

        // Solve the linear system.
        this->solve();

        // Post process of the primal and dual variables.
        this->computePrimalAndDual();

        // Compute the error.
637
        M_fixedPointResidual = ( *(this->M_primal) - *M_primalPreviousIteration ).norm2() / this->M_primal->norm2();
638
639

        // The leader process prints the iteration data.
640
        this->M_displayer.leaderPrint( "Fixed point scheme           \n" );
641

642
643
        // Print the maximum number of iterations
        this->M_displayer.leaderPrint( "Maximum number of iterations ",
Michel Kern's avatar
Michel Kern committed
644
                                       fixedPointMaxIteration(), "\n" );
645

646
647
        // Print the actual iterations
        this->M_displayer.leaderPrint( "Iteration number             ",
Michel Kern's avatar
Michel Kern committed
648
                                       fixedPointNumIteration(), "\n" );
649

650
        // Print the tollerance
651
        this->M_displayer.leaderPrint( "Tolerance                   ",
Michel Kern's avatar
Michel Kern committed
652
                                       fixedPointTolerance(), "\n" );
653

654
655
        // Print the error reached
        this->M_displayer.leaderPrint( "Error                        ",
Michel Kern's avatar
Michel Kern committed
656
                                       fixedPointResidual(), "\n" );
657

658
    }
659

Michel Kern's avatar
Michel Kern committed
660
661
    // Check if the fixed point method reach the tolerance.
    ASSERT( fixedPointMaxIteration() > fixedPointNumIteration(), "Attention the fixed point scheme did not reach convergence." );
662
663
664

} // fixedPointScheme

665
666
667
668
// Set the first value for the fixed point method.
template<typename Mesh, typename SolverType>
void
DarcySolverNonLinear<Mesh, SolverType>::
fumagalli's avatar
fumagalli committed
669
setPrimalZeroIteration ( const function_Type& primalZeroIteration )
670
671
672
673
674
675
676
{
    // Set the function for the first iteration.
    M_primalZeroIteration = primalZeroIteration;

    // Interpolate the primal variable for the first iteration.
    this->M_primal_FESpace.interpolate( M_primalZeroIteration,
                                        *(this->M_primal),
677
                                        this->M_data.dataTime()->initialTime() );
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702

} // SetZeroIterationPrimal

// ===================================================
// Protected Methods
// ===================================================

// Update all the variables of the problem.
template<typename Mesh, typename SolverType>
void
DarcySolverNonLinear<Mesh, SolverType>::
updateVariables ()
{
    // Reset the primal vector at the previous iteration step.
    M_primalPreviousIteration.reset( new vector_Type( this->M_primal_FESpace.map() ) );

    // Update the primal vector at the previous iteration step.
    *M_primalPreviousIteration = *(this->M_primal);

    // Call the method of the DarcySolver to update all the variables defined in it.
    DarcySolver<Mesh, SolverType>::updateVariables();

} // updateVariables


703
704
705
} // namespace LifeV


706
#endif //_DARCYSOLVERNONLINEAR_H_