OseenSolver.hpp 73 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
//@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 LifeV.

    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 LifeV.  If not, see <http://www.gnu.org/licenses/>.

*******************************************************************************
*/
//@HEADER

/*!
    @file
    @brief This file contains an Oseen equation solver class.

    @author Gilles Fourestey <gilles.fourestey@imag.fr>
            Simone Deparis   <simone.deparis@epfl.ch>
    @contributor Zhen Wang <zhen.wang@emory.edu>

    @date 01-06-2007

    This file contains an Oseen equation solver class.
    The resulting linear systems are solved by GMRES on the full
    matrix ( u and p coupled ).

 */


#ifndef OSEENSOLVER_H
#define OSEENSOLVER_H 1

#include <lifev/core/algorithm/SolverAztecOO.hpp>
#include <lifev/core/algorithm/Preconditioner.hpp>
#include <lifev/core/algorithm/PreconditionerIfpack.hpp>
#include <lifev/core/algorithm/PreconditionerAztecOO.hpp>
#include <lifev/core/array/MapEpetra.hpp>

#include <lifev/core/array/MatrixElemental.hpp>
#include <lifev/core/array/VectorElemental.hpp>
#include <lifev/core/array/MatrixEpetra.hpp>
#include <lifev/core/array/VectorEpetra.hpp>

#include <lifev/core/util/LifeChrono.hpp>

#include <lifev/core/fem/Assembly.hpp>
#include <lifev/core/fem/BCManage.hpp>
#include <lifev/core/fem/AssemblyElemental.hpp>
#include <lifev/core/fem/SobolevNorms.hpp>
#include <lifev/core/fem/GeometricMap.hpp>
#include <lifev/core/fem/PostProcessingBoundary.hpp>
#include <lifev/core/fem/FESpace.hpp>

#include <lifev/navier_stokes/solver/StabilizationIP.hpp>
#include <lifev/navier_stokes/solver/OseenData.hpp>

#include <boost/shared_ptr.hpp>

#include <list>

namespace LifeV
{
//! @class Oseen
/*!
    @brief This class contains an Oseen equation solver.

    @author Gilles Fourestey <gilles.fourestey@imag.fr>
            Simone Deparis   <simone.deparis@epfl.ch>
    @contributor Zhen Wang <zhen.wang@emory.edu>

 */

template< typename MeshType, typename SolverType = LifeV::SolverAztecOO >
class OseenSolver
{

public:

    //! @name Public Types
    //@{

    typedef MeshType                                    mesh_Type;
    typedef SolverType                                  linearSolver_Type;
    typedef boost::shared_ptr<linearSolver_Type>        linearSolverPtr_Type;
    typedef OseenData                                   data_Type;
    typedef boost::shared_ptr< data_Type >              dataPtr_Type;

Alessio Fumagalli's avatar
Alessio Fumagalli committed
102
103
    typedef boost::function < Real ( const Real& t, const Real& x, const Real& y,
                                     const Real& z, const ID& i ) > function_Type;
104

Alessio Fumagalli's avatar
Alessio Fumagalli committed
105
106
    typedef boost::function < Real ( const Real& t, const Real& x, const Real& y,
                                     const Real& z, const ID& i ) > source_Type;
107
108
109
110
111
112
113

    typedef BCHandler                                   bcHandler_Type;
    typedef boost::shared_ptr<bcHandler_Type>           bcHandlerPtr_Type;

#ifdef HAVE_NS_PREC
    typedef MatrixEpetraStructured<Real>                matrix_Type;
#else
114
    typedef typename linearSolver_Type::matrix_type     matrix_Type;
115
116
117
118
119
#endif
    typedef boost::shared_ptr<matrix_Type>              matrixPtr_Type;
    typedef typename linearSolver_Type::vector_type     vector_Type;
    typedef boost::shared_ptr<vector_Type>              vectorPtr_Type;

120
121
122
123
124
    typedef vector_Type                                 solution_Type;
    typedef boost::shared_ptr<solution_Type>            solutionPtr_Type;

    typedef typename linearSolver_Type::prec_raw_type   preconditioner_Type;
    typedef typename linearSolver_Type::prec_type       preconditionerPtr_Type;
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142

    //@}

    //! @name Constructors & Destructor
    //@{

    //! Empty constructor
    OseenSolver();

    //! Constructor
    /*!
        @param dataType OseenData class
        @param velocityFESpace Velocity FE space
        @param pressureFESpace Pressure FE space
        @param communicator MPI communicator
        @param lagrangeMultiplier Lagrange multiplier
     */

Alessio Fumagalli's avatar
Alessio Fumagalli committed
143
144
145
146
147
    OseenSolver ( boost::shared_ptr<data_Type>    dataType,
                  FESpace<mesh_Type, MapEpetra>&  velocityFESpace,
                  FESpace<mesh_Type, MapEpetra>&  pressureFESpace,
                  boost::shared_ptr<Epetra_Comm>& communicator,
                  const Int                       lagrangeMultiplier = 0 );
148
149
150
151
152
153
154
155
156
157

    //! Constructor
    /*!
        @param dataType OseenData class
        @param velocityFESpace Velocity FE space
        @param pressureFESpace Pressure FE space
        @param communicator MPI communicator
        @param monolithicMap MapEpetra class
        @param offset
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
158
159
160
161
162
163
    OseenSolver ( boost::shared_ptr<data_Type>    dataType,
                  FESpace<mesh_Type, MapEpetra>&  velocityFESpace,
                  FESpace<mesh_Type, MapEpetra>&  pressureFESpace,
                  boost::shared_ptr<Epetra_Comm>& communicator,
                  const MapEpetra                 monolithicMap,
                  const UInt                      offset = 0 );
164
165
166
167
168
169
170
171
172

    //! Constructor
    /*!
        @param dataType OseenData class
        @param velocityFESpace Velocity FE space
        @param pressureFESpace Pressure FE space
        @param lagrangeMultipliers (lagrange multipliers for the flux problem with rufaec flag)
        @param communicator MPI communicator
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
173
174
175
176
177
    OseenSolver ( boost::shared_ptr<data_Type>    dataType,
                  FESpace<mesh_Type, MapEpetra>&  velocityFESpace,
                  FESpace<mesh_Type, MapEpetra>&  pressureFESpace,
                  const std::vector<Int>&         lagrangeMultipliers,
                  boost::shared_ptr<Epetra_Comm>& communicator );
178
179
180
181
182
183
184
185
186
187
188
189
190

    //! virtual destructor
    virtual ~OseenSolver();

    //@}

    //! @name Methods
    //@{

    //! Set up data from GetPot
    /*!
        @param dataFile GetPot object
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
191
    virtual void setUp ( const GetPot& dataFile );
192
193
194
195
196
197

    //! Initialize with velocityFunction and pressureFunction
    /*!
        @param velocityFunction
        @param pressureFunction
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
198
    void initialize ( const function_Type& velocityFunction, const function_Type& pressureFunction );
199
200
201
202
203
204

    //! Initialize with velocityInitialGuess and pressureInitialGuess
    /*!
        @param velocityInitialGuess
        @param pressureInitialGuess
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
205
    void initialize ( const vector_Type& velocityInitialGuess, const vector_Type& pressureInitialGuess );
206
207
208
209
210

    //! Initialize with velocityAndPressure
    /*!
        @param velocityAndPressure
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
211
    void initialize ( const vector_Type& velocityAndPressure );
212
213
214
215
216
217
218
219
220
221

    //! Build linear system.
    virtual void buildSystem();

    //! Update system
    /*!
        @param alpha
        @param betaVector
        @param sourceVector
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
222
223
224
    virtual void updateSystem ( const Real         alpha,
                                const vector_Type& betaVector,
                                const vector_Type& sourceVector );
225
226
227
228
229
230
231
232
233

    //! Update system
    /*!
        @param alpha
        @param betaVector
        @param sourceVector
        @param matrix
        @param un
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
234
235
236
237
238
    virtual void updateSystem ( const Real         alpha,
                                const vector_Type& betaVector,
                                const vector_Type& sourceVector,
                                matrixPtr_Type     matrix,
                                const vector_Type&     un );
239
240
241
242
243

    //! Update stabilization term
    /*!
        @param matrixFull
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
244
    void updateStabilization ( matrix_Type& matrixFull );
245
246
247
248
249

    //! Update the right hand side
    /*!
        @param rightHandSide right hand side
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
250
    virtual void updateRightHandSide ( const vector_Type& rightHandSide )
251
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
252
        M_rightHandSideNoBC.reset (new vector_Type (rightHandSide) );
253
254
255
256
257
258
259
        M_rightHandSideNoBC->globalAssemble();
    }

    //! Update the source term
    /*!
     @param sourceTerm
    */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
260
    void updateSourceTerm (const source_Type& source );
261
262
263
264
265

    //! Update convective term, boundary condition and solve the linearized ns system
    /*!
        @param bcHandler BC handler
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
266
    virtual void iterate ( bcHandler_Type& bcHandler );
267
268
269
270
271
272

    //! Reduce the local solution in global vectors
    /*!
        @param velocity
        @param pressure
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
273
    void reduceSolution ( Vector& velocity, Vector& pressure );
274
275
276
277
278

    //! Reduce the residual
    /*!
        @param residual
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
279
    void reduceResidual ( Vector& residual );
280
281
282
283
284

    //! Set a block preconditioner
    /*!
        @blockPrecconditioner Block preconditioner
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
285
    void setBlockPreconditioner ( matrixPtr_Type blockPreconditioner );
286
287
288
289
290

    //! Update and return the coefficient matrix
    /*!
        @param matrixFull The coefficient matrix
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
291
    void getFluidMatrix ( matrix_Type& matrixFull );
292
293
294
295
296
297
298
299
300

    //! Set up post processing
    void setupPostProc( );

    //! Compute area on a boundary face with given flag
    /*!
        @param  flag
        @return area
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
301
    Real area ( const markerID_Type& flag );
302

303
304
305
306
307
    //! Compute the outgoing normal of a boundary face with given flag
    /*!
        @param  flag
        @return boundary normal
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
308
    Vector normal ( const markerID_Type& flag );
309
310
311
312
313
314

    //! Compute the geometric center of a boundary face with given flag
    /*!
        @param  flag
        @return geometric center
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
315
    Vector geometricCenter ( const markerID_Type& flag );
316

317
318
319
320
321
322
    //! Compute flux on a boundary face with given flag and a given solution
    /*!
        @param  flag
        @param  solution
        @return flux
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
323
    Real flux ( const markerID_Type& flag, const vector_Type& solution );
324
325
326
327
328
329

    //! Compute flux on a boundary face with given flag
    /*!
        @param flag
        @return flux
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
330
    Real flux ( const markerID_Type& flag );
331

332
    //! Compute the kinetic normal stress (i.e., the normal stress due to the kinetic energy) on a boundary face with a given flag and a given solution
333
    /*!
334
335
336
337
338
339
340
341
     *  @see \cite BlancoMalossi2012 \cite Malossi-Thesis
     *
     *  This method computes the following quantity:
     *
     *  \f[
     *  \mathcal{K} = \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left({\mathbf{u}}_\textrm{F} \mathbf{\cdot} {\mathbf{n}}_\textrm{F}\right)^2  \textrm{d} \Gamma
     *  \f]
     *
Cristiano Malossi's avatar
Cristiano Malossi committed
342
343
     *  @param flag boundary flag
     *  @param solution problem solution
344
     *  @return kinetic normal stress
345
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
346
    Real kineticNormalStress ( const markerID_Type& flag, const vector_Type& solution );
347

348
    //! Compute the kinetic normal stress (i.e., the normal stress due to the kinetic energy) on a boundary face with a given flag
349
    /*!
350
351
352
353
354
355
356
357
     *  @see \cite BlancoMalossi2012 \cite Malossi-Thesis
     *
     *  This method computes the following quantity:
     *
     *  \f[
     *  \mathcal{K} = \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left({\mathbf{u}}_\textrm{F} \mathbf{\cdot} {\mathbf{n}}_\textrm{F}\right)^2  \textrm{d} \Gamma
     *  \f]
     *
Cristiano Malossi's avatar
Cristiano Malossi committed
358
     *  @param flag boundary flag
359
     *  @return kinetic normal stress
360
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
361
    Real kineticNormalStress ( const markerID_Type& flag );
362

363
364
365
366
367
368
    //! Compute average pressure on a boundary face with given flag and a given solution
    /*!
        @param  flag
        @param  solution
        @return average pressure
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
369
    Real pressure ( const markerID_Type& flag, const vector_Type& solution );
370
371
372
373
374
375

    //! Compute average pressure on a boundary face with given flag
    /*!
        @param flag
        @return average pressure
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
376
    Real pressure ( const markerID_Type& flag );
377

378
    //! Compute the mean normal stress on a boundary face with a given flag and a given solution
379
    /*!
380
381
382
383
384
385
386
387
388
389
390
391
392
393
     *  @see \cite BlancoMalossi2012 \cite Malossi-Thesis
     *
     *  The mean normal stress is defined as the average of the normal component of the traction vector.
     *
     *  \f[
     *  \mathcal{S}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}
     *  \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma
     *  \f]
     *
     *  @param flag Flag of the boundary face
     *  @param bcHandler BChandler containing the boundary conditions of the problem.
     *  @param solution Vector containing the solution of the problem
     *                   (and also the Lagrange multipliers at the end).
     *  @return mean normal stress
394
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
395
    Real meanNormalStress ( const markerID_Type& flag, bcHandler_Type& bcHandler, const vector_Type& solution );
396

397
    //! Compute the mean normal stress on a boundary face with a given flag
398
    /*!
399
400
401
402
403
404
405
406
407
     *  @see \cite BlancoMalossi2012 \cite Malossi-Thesis
     *
     *  The mean normal stress is defined as the average of the normal component of the traction vector.
     *
     *  \f[
     *  \mathcal{S}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}
     *  \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma
     *  \f]
     *
408
409
410
411
     *  TODO The current version returns the exact mean normal stress if a flow rate boundary condition is imposed on the chosen boundary face.
     *  On the contrary, if other boundary conditions are applied, the mean normal stress is approximated with the mean pressure, which is a
     *  reasonable approximation for several applications.
     *
412
413
414
     *  @param flag Flag of the boundary face
     *  @param bcHandler BChandler containing the boundary conditions of the problem.
     *  @return mean normal stress
415
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
416
    Real meanNormalStress ( const markerID_Type& flag, bcHandler_Type& bcHandler );
417

418
    //! Compute the mean total normal stress on a boundary face with a given flag and a given solution
419
    /*!
420
421
422
423
424
425
426
427
428
429
     *  @see \cite BlancoMalossi2012 \cite Malossi-Thesis
     *
     *  The mean total normal stress is defined as the average of the normal component of the traction vector minus the kinetic contribution.
     *
     *  \f[
     *  \mathcal{T}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}
     *  \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma
     *  - \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left(\mathbf{u}_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right)^2 \: \textrm{d} \Gamma
     *  \f]
     *
430
431
432
433
     *  TODO The current version returns the exact mean normal stress if a flow rate boundary condition is imposed on the chosen boundary face.
     *  On the contrary, if other boundary conditions are applied, the mean normal stress is approximated with the mean pressure, which is a
     *  reasonable approximation for several applications.
     *
434
435
436
437
438
     *  @param flag Flag of the boundary face
     *  @param bcHandler BChandler containing the boundary conditions of the problem.
     *  @param solution Vector containing the solution of the problem
     *                   (and also the Lagrange multipliers at the end).
     *  @return mean total normal stress
439
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
440
    Real meanTotalNormalStress ( const markerID_Type& flag, bcHandler_Type& bcHandler, const vector_Type& solution );
441

442
    //! Compute the mean total normal stress on a boundary face with a given flag
443
    /*!
444
445
446
447
448
449
450
451
452
453
454
455
456
     *  @see \cite BlancoMalossi2012 \cite Malossi-Thesis
     *
     *  The mean total normal stress is defined as the average of the normal component of the traction vector minus the kinetic contribution.
     *
     *  \f[
     *  \mathcal{T}^{3-D}_j = \frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}
     *  \left(\sigma_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right) \mathbf{\cdot} \mathbf{n}_\textrm{F} \: \textrm{d} \Gamma
     *  - \frac{1}{2}\rho_\textrm{F}\frac{1}{\left|\Gamma^t_{\textrm{F},j}\right|}\displaystyle\int_{\Gamma^t_{\textrm{F},j}}\left(\mathbf{u}_\textrm{F} \mathbf{\cdot} \mathbf{n}_\textrm{F}\right)^2 \: \textrm{d} \Gamma
     *  \f]
     *
     *  @param flag Flag of the boundary face
     *  @param bcHandler BChandler containing the boundary conditions of the problem.
     *  @return mean total normal stress
457
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
458
    Real meanTotalNormalStress ( const markerID_Type& flag, bcHandler_Type& bcHandler );
459

460
461
462
463
464
465
466
    //! Get the Lagrange multiplier related to a flux imposed on a given part of the boundary
    /*!
        @param flag      Flag of the boundary face associated with the flux
                         and the Lagrange multiplier we want.
        @param bcHandler BChandler containing the boundary conditions of the problem.
        @return          Lagrange multiplier
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
467
    Real lagrangeMultiplier ( const markerID_Type& flag, bcHandler_Type& bcHandler );
468
469
470
471
472
473
474
475
476
477

    //! Get the Lagrange multiplier related to a flux imposed on a given part of the boundary
    /*!
        @param flag      Flag of the boundary face associated
                         with the flux and the Lagrange multiplier we want.
        @param bcHandler BChandler containing the boundary conditions of the problem.
        @param solution  Vector containing the solution of the problem
                         (and also the Lagrange multipliers at the end).
        @return          Lagrange multiplier
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
478
479
480
    Real lagrangeMultiplier ( const markerID_Type&  flag,
                              bcHandler_Type& bcHandler,
                              const vector_Type& solution );
481
482
483
484
485

    //! Reset the preconditioner.
    /*!
        @param reset Reset preconditioner.
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
486
    void resetPreconditioner ( bool reset = true )
487
488
    {
        if ( reset )
Alessio Fumagalli's avatar
Alessio Fumagalli committed
489
        {
490
            M_linearSolver->resetPreconditioner();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
491
        }
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
    }

    //! Reset stabilization matrix at the same time as the preconditioner
    void resetStabilization()
    {
        M_resetStabilization = true;
    }

    //! Update
    void updateUn()
    {
        *M_un = *M_solution;
    }

    //! Update for the monolithic
Alessio Fumagalli's avatar
Alessio Fumagalli committed
507
    void updateUn ( const vector_Type& solution )
508
509
510
511
512
513
514
515
    {
        *M_un = solution;
    }

    //! Display general information about the content of the class
    /*!
        @param output specify the output format (std::cout by default)
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
516
    void showMe ( std::ostream& output = std::cout ) const;
517
518
519
520
521
522
523
524
525
526

    //@}

    //! @name Set Methods
    //@{

    //! Set
    /*!
        @param recomputeMatrix
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
527
    void setRecomputeMatrix ( const bool& recomputeMatrix )
528
529
530
531
532
533
534
535
    {
        M_recomputeMatrix = recomputeMatrix;
    }

    //! set the source term functor
    /*!
        @param source
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
536
    void setSourceTerm ( source_Type source )
537
538
539
540
541
542
543
544
545
    {
        M_source = source;
    }

    //! Set the tolerance and the maximum number of iterations of the linear solver
    /*!
        @param tolerance Tolerance
        @param maxIteration maximum number of iterations
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
546
    void setTolMaxIteration ( const Real& tolerance, const Int& maxIteration = -1 );
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
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
703
704
705
706
707
708
709
710
711
712

    //@}

    //! @name Get Methods
    //@{

    //! Return the data container of the fluid
    /*!
        @return data container of the fluid
     */
    const dataPtr_Type& data() const
    {
        return M_oseenData;
    }

    //! Return the density of the fluid
    /*!
        @return Density of the fluid
     */
    const Real& density() const
    {
        return M_oseenData->density();
    }

    //! Return the viscosity of the fluid
    /*!
        @return Viscosity of the fluid
     */
    const Real& viscosity() const
    {
        return M_oseenData->viscosity();
    }

    //! Return the local solution vector
    /*!
        @return vectorPtr_Type Solution vector
     */
    const vectorPtr_Type& solution() const
    {
        return M_solution;
    }

    //! Return the local residual vector
    /*!
        @return Residual vector
     */
    const vector_Type& residual() const
    {
        return *M_residual;
    }

    //! Return velocity FE space
    /*!
        @return velocity FE space
     */
    FESpace<mesh_Type, MapEpetra>& velocityFESpace()
    {
        return M_velocityFESpace;
    }

    const FESpace<mesh_Type, MapEpetra>& velocityFESpace() const
    {
        return M_velocityFESpace;
    }

    //! Return pressure FE space
    /*!
        @return pressure FE space
     */
    FESpace<mesh_Type, MapEpetra>& pressureFESpace()
    {
        return M_pressureFESpace;
    }

    const FESpace<mesh_Type, MapEpetra>& pressureFESpace() const
    {
        return M_pressureFESpace;
    }

    //! Get the source term
    /*!
        @return Source term
     */
    const source_Type& sourceTerm() const
    {
        return M_source;
    }

    //! Returns the post processing structure
    /*!
        @return Post processing
    */
    PostProcessingBoundary<mesh_Type>& postProcessing()
    {
        return *M_postProcessing;
    }

    const PostProcessingBoundary<mesh_Type>& postProcessing() const
    {
        return *M_postProcessing;
    }

    //! Return MapEpetra.
    /*!
        @return MapEpetra
     */
    const MapEpetra& getMap() const
    {
        return M_localMap;
    }

    //! Return Epetra communicator
    /*!
        @return Epetra communicator
     */
    const boost::shared_ptr<Epetra_Comm>& comm() const
    {
        return M_Displayer.comm();
    }

    //! Return displayer
    /*!
        @return
     */
    const Displayer& getDisplayer() const
    {
        return M_Displayer;
    }

    //! Return
    /*!
        @return recomputeMatrix
     */
    const bool& recomputeMatrix() const
    {
        return M_recomputeMatrix;
    }

    //! Return matrix without boundary conditions
    /*!
        @return Matrix without boundary conditions
     */
    matrix_Type& matrixNoBC()
    {
        return *M_matrixNoBC;
    }

    const matrix_Type& matrixNoBC() const
    {
        return *M_matrixNoBC;
    }

    //! Return mass matrix
    /*!
        @return Mass matrix
     */
    matrix_Type& matrixMass()
    {
        return *M_velocityMatrixMass;
    }

    const matrix_Type& matrixMass() const
    {
        return *M_velocityMatrixMass;
    }

713
714
715
716
717
    const matrixPtr_Type matrixMassPtr() const
    {
        return M_velocityMatrixMass;
    }

718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
    //@}

    //@{ unused methods

    //! Set up post processing structures
    void postProcessingSetArea();

    //! Set up post processing
    void postProcessingSetNormal();

    //! Set up post processing
    void postProcessingSetPhi();

    //! Return a bool value if using diagonal block preconditioner
    bool getIsDiagonalBlockPreconditioner()
    {
        return M_isDiagonalBlockPreconditioner;
    }

    const bool& getIsDiagonalBlockPreconditioner() const
    {
        return M_isDiagonalBlockPreconditioner;
    }

    //@}

    //! Return a shared pointer to the preconditioner (of type derived from EpetraPreconditioner)
Alessio Fumagalli's avatar
Alessio Fumagalli committed
745
746
747
748
    preconditionerPtr_Type& preconditioner()
    {
        return M_linearSolver.preconditioner();
    }
749
750
751
752
753
754
755

protected:

    //! @name Constructor
    //@{

    //! Empty copy constructor
Alessio Fumagalli's avatar
Alessio Fumagalli committed
756
    OseenSolver ( const OseenSolver& oseen);
757
758
759
760
761
762
763
764
765
766
767

    //@}

    //! @name Private Methods
    //@{

    //! Removes mean of component of vector x
    /*!
        @param x
        @return
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
768
    Real removeMean ( vector_Type& x );
769
770
771
772
773
774
775

    //! Apply boundary conditions.
    /*!
        @param matrix
        @param rightHandSide
        @param bcHandler
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
776
777
778
    void applyBoundaryConditions ( matrix_Type&        matrix,
                                   vector_Type&        rightHandSide,
                                   bcHandler_Type& bcHandler );
779
780
781
782
783

    //! Echo message.
    /*!
        @param message
     */
Alessio Fumagalli's avatar
Alessio Fumagalli committed
784
    void echo ( std::string message );
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823

    //! Return the dim of velocity FE space
    const UInt& dimVelocity() const
    {
        return M_velocityFESpace.dim();
    }

    //! Return the dim of pressure FE space
    const UInt& dimPressure() const
    {
        return M_pressureFESpace.dim();
    }

    //@}

    //private members

    //! data for Navier-Stokes solvers
    dataPtr_Type                   M_oseenData;

    // FE spaces
    FESpace<mesh_Type, MapEpetra>& M_velocityFESpace;
    FESpace<mesh_Type, MapEpetra>& M_pressureFESpace;

    //! MPI communicator
    Displayer                      M_Displayer;

    MapEpetra                      M_localMap;

    //! mass matrix
    matrixPtr_Type                 M_velocityMatrixMass;

    //! mass matrix
    matrixPtr_Type                 M_pressureMatrixMass;

    //! Stokes matrix: nu*stiff
    matrixPtr_Type                 M_matrixStokes;

    //! matrix to be solved
Alessio Fumagalli's avatar
Alessio Fumagalli committed
824
    //    matrixPtr_Type               M_matrixFull;
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898

    //! matrix without boundary conditions
    matrixPtr_Type                 M_matrixNoBC;

    //! stabilization matrix
    matrixPtr_Type                 M_matrixStabilization;

    //! source term for Navier-Stokes equations
    source_Type                    M_source;

    //! Right hand side for the velocity component
    vectorPtr_Type                 M_rightHandSideNoBC;

    //! Global solution
    vectorPtr_Type                 M_solution;

    //! residual
    vectorPtr_Type                    M_residual;

    linearSolverPtr_Type              M_linearSolver;

    bool                           M_steady;

    //! Postprocessing class
    boost::shared_ptr<PostProcessingBoundary<mesh_Type> > M_postProcessing;

    //! Stabilization
    bool                           M_stabilization;
    bool                           M_reuseStabilization;
    bool                           M_resetStabilization;
    Int                            M_iterReuseStabilization;

    details::StabilizationIP<mesh_Type, DOF> M_ipStabilization;
    Real                           M_gammaBeta;
    Real                           M_gammaDiv;
    Real                           M_gammaPress;

    const function_Type*                M_betaFunction;

    bool                           M_divBetaUv;

    bool                           M_stiffStrain;

    //
    Real                           M_diagonalize;

    UInt                           M_count;

    bool                           M_recomputeMatrix;

    bool                           M_isDiagonalBlockPreconditioner;

    //! Elementary matrices and vectors
    MatrixElemental                        M_elementMatrixStiff;      // velocity Stokes
    MatrixElemental                        M_elementMatrixMass;       // velocity mass
    MatrixElemental                        M_elementMatrixPreconditioner;          // (p,q) bloc for preconditioners
    MatrixElemental                        M_elementMatrixDivergence;
    MatrixElemental                        M_elementMatrixGradient;
    VectorElemental                        M_elementRightHandSide;           // Elementary right hand side
    matrixPtr_Type                         M_blockPreconditioner;
    VectorElemental                        M_wLoc;
    VectorElemental                        M_uLoc;
    boost::shared_ptr<vector_Type> M_un;

}; // class OseenSolver



// ===================================================
// Constructors & Destructor
// ===================================================

template<typename MeshType, typename SolverType>
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
OseenSolver ( boost::shared_ptr<data_Type>    dataType,
              FESpace<mesh_Type, MapEpetra>&  velocityFESpace,
              FESpace<mesh_Type, MapEpetra>&  pressureFESpace,
              boost::shared_ptr<Epetra_Comm>& communicator,
              const Int                       lagrangeMultiplier ) :
    M_oseenData       ( dataType ),
    M_velocityFESpace        ( velocityFESpace ),
    M_pressureFESpace        ( pressureFESpace ),
    M_Displayer              ( communicator ),
    M_localMap               ( M_velocityFESpace.map() + M_pressureFESpace.map() + lagrangeMultiplier),
    M_velocityMatrixMass     ( ),
    M_pressureMatrixMass     ( ),
    M_matrixStokes           ( ),
    M_matrixNoBC             ( ),
    M_matrixStabilization    ( ),
    M_rightHandSideNoBC      ( ),
    M_solution               ( new vector_Type ( M_localMap ) ),
    M_residual               ( new vector_Type (M_localMap ) ),
    M_linearSolver           ( new linearSolver_Type (communicator) ),
    M_steady                 ( ),
    M_postProcessing         ( new PostProcessingBoundary<mesh_Type> ( M_velocityFESpace.mesh(),
                                                                       &M_velocityFESpace.feBd(),
                                                                       &M_velocityFESpace.dof(),
                                                                       &M_pressureFESpace.feBd(),
                                                                       &M_pressureFESpace.dof(),
                                                                       M_localMap ) ),
    M_stabilization          ( false ),
    M_reuseStabilization     ( false ),
    M_resetStabilization     ( false ),
    M_iterReuseStabilization ( -1 ),
    //        M_ipStabilization        ( M_velocityFESpace.mesh(),
    //                                   M_velocityFESpace.dof(),
    //                                   M_velocityFESpace.refFE(),
    //                                   M_velocityFESpace.feBd(),
    //                                   M_velocityFESpace.qr(),
    //                                   0., 0., 0.,
    //                                   M_oseenData->viscosity() ),
    M_betaFunction           ( 0 ),
    M_divBetaUv              ( false ),
    M_stiffStrain            ( false ),
    M_diagonalize            ( false ),
    M_count                  ( 0 ),
    M_recomputeMatrix        ( false ),
    M_elementMatrixStiff     ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim(), velocityFESpace.fieldDim() ),
    M_elementMatrixMass      ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim(), velocityFESpace.fieldDim() ),
    M_elementMatrixPreconditioner ( M_pressureFESpace.fe().nbFEDof(), 1, 1 ),
    M_elementMatrixDivergence ( M_pressureFESpace.fe().nbFEDof(), 1, 0,
                                M_velocityFESpace.fe().nbFEDof(), 0, velocityFESpace.fieldDim() ),
    M_elementMatrixGradient  ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim(), 0,
                               M_pressureFESpace.fe().nbFEDof(), 0, 1 ),
    M_elementRightHandSide   ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim() ),
    M_blockPreconditioner    ( ),
    M_wLoc                   ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim() ),
    M_uLoc                   ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim() ),
    M_un                     ( new vector_Type (M_localMap) )
954
955
956
{
    // if(M_stabilization = ( &M_velocityFESpace.refFE() == &M_pressureFESpace.refFE() ))
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
957
958
        M_ipStabilization.setFeSpaceVelocity (M_velocityFESpace);
        M_ipStabilization.setViscosity (M_oseenData->viscosity() );
959
960
961
962
963
    }
}

template<typename MeshType, typename SolverType>
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
OseenSolver ( boost::shared_ptr<data_Type>    dataType,
              FESpace<mesh_Type, MapEpetra>&  velocityFESpace,
              FESpace<mesh_Type, MapEpetra>&  pressureFESpace,
              boost::shared_ptr<Epetra_Comm>& communicator,
              MapEpetra                       monolithicMap,
              UInt                            /*offset*/ ) :
    M_oseenData              ( dataType ),
    M_velocityFESpace        ( velocityFESpace ),
    M_pressureFESpace        ( pressureFESpace ),
    M_Displayer              ( communicator ),
    M_localMap               ( monolithicMap ),
    M_velocityMatrixMass     ( ),
    M_matrixStokes           ( ),
    M_matrixNoBC             ( ),
    M_matrixStabilization    ( ),
    M_rightHandSideNoBC      ( ),
    M_solution               ( ),
    M_residual               (  ),
    M_linearSolver           ( ),
    M_postProcessing         ( new PostProcessingBoundary<mesh_Type> (M_velocityFESpace.mesh(),
                                                                      &M_velocityFESpace.feBd(),
                                                                      &M_velocityFESpace.dof(),
                                                                      &M_pressureFESpace.feBd(),
                                                                      &M_pressureFESpace.dof(),
                                                                      M_localMap ) ),
    M_stabilization          ( false ),
    M_reuseStabilization     ( false ),
    M_resetStabilization     ( false ),
    M_iterReuseStabilization ( -1 ),
    M_betaFunction           ( 0 ),
    M_divBetaUv              ( false ),
    M_stiffStrain            ( false ),
    M_diagonalize            ( false ),
    M_count                  ( 0 ),
    M_recomputeMatrix        ( false ),
    M_elementMatrixStiff     ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim(), M_velocityFESpace.fieldDim() ),
    M_elementMatrixMass      ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim(), M_velocityFESpace.fieldDim() ),
    M_elementMatrixPreconditioner                 ( M_pressureFESpace.fe().nbFEDof(), 1, 1 ),
    M_elementMatrixDivergence ( M_pressureFESpace.fe().nbFEDof(), 1, 0,
                                M_velocityFESpace.fe().nbFEDof(), 0, M_velocityFESpace.fieldDim() ),
    M_elementMatrixGradient  ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim(), 0,
                               M_pressureFESpace.fe().nbFEDof(), 0, 1 ),
    M_elementRightHandSide   ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim() ),
    M_blockPreconditioner    ( ),
    M_wLoc                   ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim() ),
    M_uLoc                   ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim() ),
    M_un                     ( /*new vector_Type(M_localMap)*/ )
1011
1012
1013
{
    // if(M_stabilization = ( &M_velocityFESpace.refFE() == &M_pressureFESpace.refFE() ))
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1014
1015
        M_ipStabilization.setFeSpaceVelocity (M_velocityFESpace);
        M_ipStabilization.setViscosity (M_oseenData->viscosity() );
1016
1017
1018
1019
1020
    }
}

template<typename MeshType, typename SolverType>
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
OseenSolver ( boost::shared_ptr<data_Type>    dataType,
              FESpace<mesh_Type, MapEpetra>&  velocityFESpace,
              FESpace<mesh_Type, MapEpetra>&  pressureFESpace,
              const std::vector<Int>&         lagrangeMultipliers,
              boost::shared_ptr<Epetra_Comm>& communicator ) :
    M_oseenData       ( dataType ),
    M_velocityFESpace        ( velocityFESpace ),
    M_pressureFESpace        ( pressureFESpace ),
    M_Displayer              ( communicator ),
    M_localMap               ( M_velocityFESpace.map() + M_pressureFESpace.map() + lagrangeMultipliers ),
    M_velocityMatrixMass     ( ),
    M_matrixStokes           ( ),
    M_matrixNoBC             ( ),
    M_matrixStabilization    ( ),
    M_rightHandSideNoBC      ( ),
    M_solution               ( new vector_Type ( M_localMap ) ),
    M_residual               (  ),
    M_linearSolver           ( new linearSolver_Type (communicator) ),
    M_postProcessing         ( new PostProcessingBoundary<mesh_Type> (M_velocityFESpace.mesh(),
                                                                      &M_velocityFESpace.feBd(),
                                                                      &M_velocityFESpace.dof(),
                                                                      &M_pressureFESpace.feBd(),
                                                                      &M_pressureFESpace.dof(),
                                                                      M_localMap ) ),
    M_stabilization          ( false ),
    M_reuseStabilization     ( false ),
    M_resetStabilization     ( false ),
    M_iterReuseStabilization ( -1 ),
    M_betaFunction           ( 0 ),
    M_divBetaUv              ( false ),
    M_stiffStrain            ( false ),
    M_diagonalize            ( false ),
    M_count                  ( 0 ),
    M_recomputeMatrix        ( false ),
    M_elementMatrixStiff     ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim(), M_velocityFESpace.fieldDim() ),
    M_elementMatrixMass      ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim(), M_velocityFESpace.fieldDim() ),
    M_elementMatrixPreconditioner ( M_pressureFESpace.fe().nbFEDof(), 1, 1 ),
    M_elementMatrixDivergence ( M_pressureFESpace.fe().nbFEDof(), 1, 0,
                                M_velocityFESpace.fe().nbFEDof(), 0, M_velocityFESpace.fieldDim() ),
    M_elementMatrixGradient  ( M_velocityFESpace.fe().nbFEDof(), M_velocityFESpace.fieldDim(), 0,
                               M_pressureFESpace.fe().nbFEDof(), 0, 1 ),
    M_elementRightHandSide   ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim() ),
    M_blockPreconditioner    ( ),
    M_wLoc                   ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim() ),
    M_uLoc                   ( M_velocityFESpace.fe().nbFEDof(), velocityFESpace.fieldDim() ),
    M_un                     ( new vector_Type (M_localMap) )
1067
1068
1069
{
    // if(M_stabilization = ( &M_velocityFESpace.refFE() == &M_pressureFESpace.refFE() ))
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1070
1071
        M_ipStabilization.setFeSpaceVelocity (M_velocityFESpace);
        M_ipStabilization.setViscosity (M_oseenData->viscosity() );
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
    }
}

template<typename MeshType, typename SolverType>
OseenSolver<MeshType, SolverType>::
~OseenSolver()
{

}


// ===================================================
// Methods
// ===================================================

template<typename MeshType, typename SolverType>
void
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1089
OseenSolver<MeshType, SolverType>::setUp ( const GetPot& dataFile )
1090
{
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1091
    if (M_linearSolver.get() )
1092
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1093
1094
        M_linearSolver->setupPreconditioner ( dataFile, "fluid/prec" );
        M_linearSolver->setDataFromGetPot ( dataFile, "fluid/solver" );
1095
1096
    }

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1097
1098
1099
    M_stabilization = dataFile ( "fluid/ipstab/use",  (&M_velocityFESpace.refFE() == &M_pressureFESpace.refFE() ) );
    M_steady        = dataFile ( "fluid/miscellaneous/steady", 0 );
    if (M_stabilization)
1100
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1101
1102
1103
1104
1105
1106
1107
        M_gammaBeta     = dataFile ( "fluid/ipstab/gammaBeta",  0. );
        M_gammaDiv      = dataFile ( "fluid/ipstab/gammaDiv",   0. );
        M_gammaPress    = dataFile ( "fluid/ipstab/gammaPress", 0. );
        M_reuseStabilization     = dataFile ( "fluid/ipstab/reuse", false );
        if (M_linearSolver.get() )
            M_iterReuseStabilization = dataFile ( "fluid/ipstab/max_iter_reuse",
                                                  static_cast<Int> ( M_linearSolver->maxNumIterations() * 8. / 10. ) );
1108
1109
    }
    // Energetic stabilization term
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1110
    M_divBetaUv   = dataFile ( "fluid/space_discretization/div_beta_u_v", false);
1111
    // Enable grad( u )^T in stress tensor
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1112
1113
1114
    M_stiffStrain = dataFile ( "fluid/space_discretization/stiff_strain", false);
    M_diagonalize = dataFile ( "fluid/space_discretization/diagonalize", 1. );
    M_isDiagonalBlockPreconditioner = dataFile ( "fluid/diagonalBlockPrec", false );
1115
1116
1117
1118
1119

    //    M_linearSolver.setAztecooPreconditioner( dataFile, "fluid/solver" );

    M_ipStabilization.setGammaBeta ( M_gammaBeta );
    M_ipStabilization.setGammaDiv  ( M_gammaDiv );
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1120
    M_ipStabilization.setGammaPress ( M_gammaPress );
1121
1122
1123
1124
1125
1126
}


template<typename MeshType, typename SolverType>
void
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1127
initialize ( const function_Type& velocityFunction, const function_Type& pressureFunction )
1128
{
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1129
1130
1131
1132
    vector_Type velocityInitialGuess ( M_velocityFESpace.map() );
    M_velocityFESpace.interpolate ( velocityFunction,
                                    velocityInitialGuess,
                                    M_oseenData->dataTime()->time() );
1133

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1134
1135
1136
1137
    vector_Type pressureInitialGuess ( M_pressureFESpace.map() );
    M_pressureFESpace.interpolate ( pressureFunction,
                                    pressureInitialGuess,
                                    M_oseenData->dataTime()->time() );
1138

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1139
    initialize ( velocityInitialGuess, pressureInitialGuess );
1140
1141
1142
1143
1144
1145
}


template<typename MeshType, typename SolverType>
void
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1146
initialize ( const vector_Type& velocityInitialGuess, const vector_Type& pressureInitialGuess )
1147
1148
1149
1150
{

    *M_solution = velocityInitialGuess;
    *M_un = velocityInitialGuess;
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1151
1152
    M_solution->add ( pressureInitialGuess, M_velocityFESpace.fieldDim() * M_velocityFESpace.dof().numTotalDof() );
    M_un->add ( pressureInitialGuess, M_velocityFESpace.fieldDim() * M_velocityFESpace.dof().numTotalDof() );
1153
1154
1155
1156
1157
1158

}

template<typename MeshType, typename SolverType>
void
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1159
initialize ( const vector_Type& velocityAndPressure )
1160
1161
1162
1163
{

    *M_un = velocityAndPressure;
    if ( M_solution.get() )
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1164
    {
1165
        *M_solution = velocityAndPressure;
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1166
    }
1167
1168
1169
1170
1171
1172
1173

}

template<typename MeshType, typename SolverType>
void
OseenSolver<MeshType, SolverType>::buildSystem()
{
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1174
1175
    M_velocityMatrixMass.reset  ( new matrix_Type ( M_localMap ) );
    M_matrixStokes.reset ( new matrix_Type ( M_localMap ) );
1176

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1177
    M_Displayer.leaderPrint ( "  F-  Computing constant matrices ...          " );
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199

    LifeChrono chrono;

    LifeChrono chronoDer;
    LifeChrono chronoStiff;
    LifeChrono chronoMass;
    LifeChrono chronoGrad;

    LifeChrono chronoStiffAssemble;
    LifeChrono chronoMassAssemble;
    LifeChrono chronoGradAssemble;
    LifeChrono chronoDivAssemble;
    LifeChrono chronoStab;
    LifeChrono chronoZero;

    // Number of velocity components
    UInt numVelocityComponent = M_velocityFESpace.fieldDim();

    // Elementary computation and matrix assembling
    // Loop on elements

    UInt velocityTotalDof   = M_velocityFESpace.dof().numTotalDof();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1200
    //    UInt pressureTotalDof = M_pressureFESpace.dof().numTotalDof();
1201
1202
1203

    if ( M_isDiagonalBlockPreconditioner == true )
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1204
        M_blockPreconditioner.reset ( new matrix_Type ( M_localMap ) );
1205
1206
1207
1208
1209
1210
1211
    }
    chrono.start();

    for ( UInt iElement = 0; iElement < M_velocityFESpace.mesh()->numElements(); iElement++ )
    {
        chronoDer.start();
        // just to provide the id number in the assem_mat_mixed
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1212
        M_pressureFESpace.fe().update ( M_velocityFESpace.mesh()->element ( iElement ) );
1213
1214
        // just to provide the id number in the assem_mat_mixed
        // M_pressureFESpace.fe().updateFirstDeriv( M_velocityFESpace.mesh()->element( iElement ) );
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1215
        M_velocityFESpace.fe().updateFirstDeriv ( M_velocityFESpace.mesh()->element ( iElement ) );
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229

        chronoDer.stop();

        chronoZero.start();
        M_elementMatrixStiff.zero();
        M_elementMatrixMass.zero();
        M_elementMatrixPreconditioner.zero();
        M_elementMatrixDivergence.zero();
        M_elementMatrixGradient.zero();
        chronoZero.stop();

        // stiffness matrix
        chronoStiff.start();
        if ( M_stiffStrain )
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1230
1231
1232
            stiff_strain ( 2.0 * M_oseenData->viscosity(),
                           M_elementMatrixStiff,
                           M_velocityFESpace.fe() );
1233
        else
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1234
1235
1236
            stiff ( M_oseenData->viscosity(),
                    M_elementMatrixStiff,
                    M_velocityFESpace.fe(), 0, 0, M_velocityFESpace.fieldDim() );
1237
1238
1239
1240
1241
1242
1243
        //stiff_div( 0.5*M_velocityFESpace.fe().diameter(), M_elementMatrixStiff, M_velocityFESpace.fe() );
        chronoStiff.stop();

        // mass matrix
        if ( !M_steady )
        {
            chronoMass.start();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1244
1245
1246
            mass ( M_oseenData->density(),
                   M_elementMatrixMass,
                   M_velocityFESpace.fe(), 0, 0, M_velocityFESpace.fieldDim() );
1247
1248
1249
1250
1251
1252
1253
1254
1255
            chronoMass.stop();
        }

        for ( UInt iComponent = 0; iComponent < numVelocityComponent; iComponent++ )
        {
            // stiffness matrix
            chronoStiffAssemble.start();
            if ( M_isDiagonalBlockPreconditioner == true )
            {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1256
1257
1258
1259
1260
1261
1262
1263
                assembleMatrix ( *M_blockPreconditioner,
                                 M_elementMatrixStiff,
                                 M_velocityFESpace.fe(),
                                 M_velocityFESpace.fe(),
                                 M_velocityFESpace.dof(),
                                 M_velocityFESpace.dof(),
                                 iComponent, iComponent,
                                 iComponent * velocityTotalDof, iComponent * velocityTotalDof);
1264
1265
1266
1267
1268
1269
1270
            }
            else
            {
                if ( M_stiffStrain ) // sigma = 0.5 * mu (grad( u ) + grad ( u )^T)
                {
                    for ( UInt jComp = 0; jComp < numVelocityComponent; jComp++ )
                    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1271
1272
1273
1274
1275
1276
1277
1278
                        assembleMatrix ( *M_matrixStokes,
                                         M_elementMatrixStiff,
                                         M_velocityFESpace.fe(),
                                         M_velocityFESpace.fe(),
                                         M_velocityFESpace.dof(),
                                         M_velocityFESpace.dof(),
                                         iComponent, jComp,
                                         iComponent * velocityTotalDof, jComp * velocityTotalDof);
1279
1280
1281
1282
1283

                    }
                }
                else // sigma = mu grad( u )
                {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1284
1285
1286
1287
1288
1289
1290
1291
                    assembleMatrix ( *M_matrixStokes,
                                     M_elementMatrixStiff,
                                     M_velocityFESpace.fe(),
                                     M_velocityFESpace.fe(),
                                     M_velocityFESpace.dof(),
                                     M_velocityFESpace.dof(),
                                     iComponent, iComponent,
                                     iComponent * velocityTotalDof, iComponent * velocityTotalDof);
1292
1293
1294
1295
1296
1297
1298
1299
                }
            }
            chronoStiffAssemble.stop();

            // mass matrix
            if ( !M_steady )
            {
                chronoMassAssemble.start();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1300
1301
1302
1303
1304
1305
1306
1307
                assembleMatrix ( *M_velocityMatrixMass,
                                 M_elementMatrixMass,
                                 M_velocityFESpace.fe(),
                                 M_velocityFESpace.fe(),
                                 M_velocityFESpace.dof(),
                                 M_velocityFESpace.dof(),
                                 iComponent, iComponent,
                                 iComponent * velocityTotalDof, iComponent * velocityTotalDof);
1308
1309
1310
1311
1312
                chronoMassAssemble.stop();
            }

            // divergence
            chronoGrad.start();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1313
1314
1315
1316
1317
            grad ( iComponent, 1.0,
                   M_elementMatrixGradient,
                   M_velocityFESpace.fe(),
                   M_pressureFESpace.fe(),
                   iComponent, 0 );
1318
1319
1320
            chronoGrad.stop();

            chronoGradAssemble.start();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1321
1322
1323
1324
1325
1326
1327
1328
            assembleMatrix ( *M_matrixStokes,
                             M_elementMatrixGradient,
                             M_velocityFESpace.fe(),
                             M_pressureFESpace.fe(),
                             M_velocityFESpace.dof(),
                             M_pressureFESpace.dof(),
                             iComponent, 0,
                             iComponent * velocityTotalDof, numVelocityComponent * velocityTotalDof );
1329
1330
1331
            chronoGradAssemble.stop();

            chronoDivAssemble.start();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1332
1333
1334
1335
1336
1337
1338
1339
1340
            assembleTransposeMatrix ( *M_matrixStokes,
                                      -1.,
                                      M_elementMatrixGradient,
                                      M_pressureFESpace.fe(),
                                      M_velocityFESpace.fe(),
                                      M_pressureFESpace.dof(),
                                      M_velocityFESpace.dof(),
                                      0 , iComponent,
                                      numVelocityComponent * velocityTotalDof, iComponent * velocityTotalDof );
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
            chronoDivAssemble.stop();
        }
    }

    //    for (UInt ii = M_velocityFESpace.fieldDim()*dimVelocity(); ii < M_velocityFESpace.fieldDim()*dimVelocity() + dimPressure(); ++ii)
    //  M_matrixStokes->set_mat_inc( ii ,ii, 0. ); not scalable!!!

    if ( M_isDiagonalBlockPreconditioner == true )
    {
        M_blockPreconditioner->globalAssemble();
        *M_matrixStokes += *M_blockPreconditioner;
    }
    comm()->Barrier();

    chrono.stop();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1356
    M_Displayer.leaderPrintMax ( "done in " , chrono.diff() );
1357

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1358
    M_Displayer.leaderPrint ( "  F-  Finalizing the matrices ...              " );
1359
1360
1361
1362
1363
1364
1365

    chrono.start();

    M_matrixStokes->globalAssemble();
    M_velocityMatrixMass->globalAssemble();

    chrono.stop();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1366
    M_Displayer.leaderPrintMax ( "done in " , chrono.diff() );
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385

    if ( false )
        std::cout << " partial times:  \n"
                  << " Der            " << chronoDer.diffCumul() << " s.\n"
                  << " Zero           " << chronoZero.diffCumul() << " s.\n"
                  << " Stiff          " << chronoStiff.diffCumul() << " s.\n"
                  << " Stiff Assemble " << chronoStiffAssemble.diffCumul() << " s.\n"
                  << " Mass           " << chronoMass.diffCumul() << " s.\n"
                  << " Mass Assemble  " << chronoMassAssemble.diffCumul() << " s.\n"
                  << " Grad           " << chronoGrad.diffCumul() << " s.\n"
                  << " Grad Assemble  " << chronoGradAssemble.diffCumul() << " s.\n"
                  << " Div Assemble   " << chronoDivAssemble.diffCumul() << " s.\n"
                  << std::endl;

}

template<typename MeshType, typename SolverType>
void
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1386
1387
1388
updateSystem ( const Real         alpha,
               const vector_Type& betaVector,
               const vector_Type& sourceVector )
1389
1390
{
    if ( M_matrixNoBC.get() )
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1391
1392
1393
    {
        M_matrixNoBC.reset ( new matrix_Type ( M_localMap, M_matrixNoBC->meanNumEntries() ) );
    }
1394
    else
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1395
1396
1397
    {
        M_matrixNoBC.reset ( new matrix_Type ( M_localMap ) );
    }
1398

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1399
    updateSystem ( alpha, betaVector, sourceVector, M_matrixNoBC, *M_un );
1400
    if ( alpha != 0. )
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1401
    {
1402
        M_matrixNoBC->globalAssemble();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1403
    }
1404
1405
1406
1407
1408
1409

}

template<typename MeshType, typename SolverType>
void
OseenSolver<MeshType, SolverType>::
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1410
1411
1412
1413
1414
updateSystem ( const Real         alpha,
               const vector_Type& betaVector,
               const vector_Type& sourceVector,
               matrixPtr_Type     matrixNoBC,
               const vector_Type&     un )
1415
1416
1417
1418
1419
1420
1421
{
    LifeChrono chrono;

    // clearing pressure mass matrix in case we need it in removeMean;
    M_pressureMatrixMass.reset( );


Alessio Fumagalli's avatar
Alessio Fumagalli committed
1422
    M_Displayer.leaderPrint ( "  F-  Updating mass term on right hand side... " );
1423
1424
1425
1426

    chrono.start();

    UInt velocityTotalDof   = M_velocityFESpace.dof().numTotalDof();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1427
    //    UInt pressureTotalDof = M_pressureFESpace.dof().numTotalDof();
1428
1429
1430

    // Right hand side for the velocity at time

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1431
    updateRightHandSide ( sourceVector );
1432
1433
1434

    chrono.stop();

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1435
    M_Displayer.leaderPrintMax ( "done in ", chrono.diff() );
1436
1437
1438
1439
1440


    //    M_updated = false;

    if ( M_recomputeMatrix )
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1441
    {
1442
        buildSystem();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1443
    }
1444

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1445
    M_Displayer.leaderPrint ( "  F-  Copying the matrices ...                 " );
1446
1447
1448
1449
1450

    chrono.start();

    if ( M_isDiagonalBlockPreconditioner == true )
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1451
1452
1453
        matrixPtr_Type tempMatrix ( M_blockPreconditioner );
        M_blockPreconditioner.reset ( new matrix_Type ( M_localMap,
                                                        M_blockPreconditioner->meanNumEntries() ) );
1454
1455
1456
1457
1458
        *M_blockPreconditioner += *tempMatrix;
    }


    chrono.stop();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1459
    M_Displayer.leaderPrintMax ( "done in " , chrono.diff() );
1460
1461
1462
1463
1464
1465
1466


    UInt numVelocityComponent = M_velocityFESpace.fieldDim();

    //! managing the convective term

    Real normInf;
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1467
    betaVector.normInf ( &normInf );
1468
1469
1470

    if ( normInf != 0. )
    {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1471
        M_Displayer.leaderPrint ( "  F-  Sharing convective term ...              " );
1472
1473
1474
1475
        chrono.start();

        // vector with repeated nodes over the processors

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1476
1477
        vector_Type betaVectorRepeated ( betaVector, Repeated );
        vector_Type unRepeated ( un, Repeated );
1478
1479
1480

        chrono.stop();

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1481
1482
        M_Displayer.leaderPrintMax ( "done in " , chrono.diff() );
        M_Displayer.leaderPrint ( "  F-  Updating the convective terms ...        " );
1483
1484
1485
1486
1487
        chrono.start();

        for ( UInt iElement = 0; iElement < M_velocityFESpace.mesh()->numElements(); ++iElement )
        {
            // just to provide the id number in the assem_mat_mixed
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1488
            M_pressureFESpace.fe().updateFirstDeriv ( M_velocityFESpace.mesh()->element ( iElement ) );
1489
            //as updateFirstDer
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1490
            M_velocityFESpace.fe().updateFirstDeriv ( M_velocityFESpace.mesh()->element ( iElement ) );
1491
1492
1493
1494
1495
1496
1497
1498

            M_elementMatrixStiff.zero();

            UInt elementID = M_velocityFESpace.fe().currentLocalId();
            // Non linear term, Semi-implicit approach
            // M_elementRightHandSide contains the velocity values in the nodes
            for ( UInt iNode = 0 ; iNode < M_velocityFESpace.fe().nbFEDof() ; iNode++ )
            {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1499
                UInt iLocal = M_velocityFESpace.fe().patternFirst ( iNode );
1500
1501
                for ( UInt iComponent = 0; iComponent < numVelocityComponent; ++iComponent )
                {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1502
                    UInt iGlobal = M_velocityFESpace.dof().localToGlobalMap ( elementID, iLocal )
1503
1504
                                   + iComponent * dimVelocity();
                    M_elementRightHandSide.vec() [ iLocal + iComponent * M_velocityFESpace.fe().nbFEDof() ]
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1505
                        = betaVectorRepeated[iGlobal];
1506
1507

                    M_uLoc.vec() [ iLocal + iComponent * M_velocityFESpace.fe().nbFEDof() ]
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1508
                        = unRepeated (iGlobal);
1509
                    M_wLoc.vec() [ iLocal + iComponent * M_velocityFESpace.fe().nbFEDof() ]
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1510
                        = unRepeated (iGlobal) - betaVectorRepeated (iGlobal);
1511
1512
1513
                }
            }

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1514
1515
1516
1517
1518
1519
1520
1521
            if (M_oseenData->conservativeFormulation() )
            {
                // ALE term: - rho div(w) u v
                mass_divw ( - M_oseenData->density(),
                            M_wLoc,
                            M_elementMatrixStiff,
                            M_velocityFESpace.fe(), 0, 0, numVelocityComponent );
            }
1522

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1523
1524
1525
1526
1527
            // ALE stab implicit: 0.5 rho div u w v
            mass_divw ( 0.5 * M_oseenData->density(),
                        M_uLoc,
                        M_elementMatrixStiff,
                        M_velocityFESpace.fe(), 0, 0, numVelocityComponent );
1528

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1529
1530
1531
1532
1533
1534
            // Stabilising term: -rho div(u^n) u v
            if ( M_divBetaUv )
                mass_divw ( -0.5 * M_oseenData->density(),
                            M_uLoc,
                            M_elementMatrixStiff,
                            M_velocityFESpace.fe(), 0, 0, numVelocityComponent );
1535
1536

            // compute local convective terms
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1537
1538
1539
1540
            advection ( M_oseenData->density(),
                        M_elementRightHandSide,
                        M_elementMatrixStiff,
                        M_velocityFESpace.fe(), 0, 0, numVelocityComponent );
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552

            // loop on components
            for ( UInt iComponent = 0; iComponent < numVelocityComponent; ++iComponent )
            {
                // compute local convective term and assembling
                // grad( 0, M_elementRightHandSide, M_elementMatrixStiff, M_velocityFESpace.fe(),
                //       M_velocityFESpace.fe(), iComponent, iComponent );
                // grad( 1, M_elementRightHandSide, M_elementMatrixStiff, M_velocityFESpace.fe(),
                //       M_velocityFESpace.fe(), iComponent, iComponent );
                // grad( 2, M_elementRightHandSide, M_elementMatrixStiff, M_velocityFESpace.fe(),
                //       M_velocityFESpace.fe(), iComponent, iComponent );

Alessio Fumagalli's avatar
Alessio Fumagalli committed
1553
1554
1555
1556
1557
1558
1559
1560
                assembleMatrix ( *matrixNoBC,
                                 M_elementMatrixStiff,
                                 M_velocityFESpace.fe(),
                                 M_velocityFESpace.fe(),
                                 M_velocityFESpace.dof(),
                                 M_velocityFESpace.dof(),
                                 iComponent, iComponent,
                                 iComponent * velocityTotalDof, iComponent * velocityTotalDof );
1561
1562
1563
1564
            }
        }

        chrono.stop();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1565
        M_Displayer.leaderPrintMax ( "done in " , chrono.diff() );
1566
1567
1568
1569

        if ( M_stabilization &&
                ( M_resetStabilization || !M_reuseStabilization || ( M_matrixStabilization.get() == 0 ) ) )
        {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1570
            M_Displayer.leaderPrint ( "  F-  Updating the stabilization terms ...     " );
1571
            chrono.start();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1572
1573
            M_matrixStabilization.reset ( new matrix_Type ( M_localMap ) );
            M_ipStabilization.apply ( *M_matrixStabilization, betaVectorRepeated, false );
1574
1575
1576
            M_matrixStabilization->globalAssemble();
            M_resetStabilization = false;
            chrono.stop();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1577
            M_Displayer.leaderPrintMax ( "done in " , chrono.diff() );
1578
1579
1580
1581
1582
1583
1584
        }

    }
    else
    {
        if ( M_stabilization )
        {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1585
            M_Displayer.leaderPrint ( "  F-  Updating the stabilization terms ...     " );
1586
1587
1588
1589
            chrono.start();

            if ( M_resetStabilization || !M_reuseStabilization || ( M_matrixStabilization.get() == 0 ) )
            {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1590
1591
                M_matrixStabilization.reset ( new matrix_Type ( M_localMap ) );
                M_ipStabilization.apply ( *M_matrixStabilization, betaVector, false );
1592
1593
1594
                M_matrixStabilization->globalAssemble();
                M_resetStabilization = false;
                chrono.stop();
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1595
                M_Displayer.leaderPrintMax ( "done in " , chrono.diff() );
1596
1597
1598
            }
            else
            {
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1599
                M_Displayer.leaderPrint ( "reusing\n" );
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
            }
        }
    }

    if ( alpha != 0. )
    {
        *matrixNoBC += (*M_velocityMatrixMass) * alpha;
        if ( M_isDiagonalBlockPreconditioner == true )
        {
            matrixNoBC->globalAssemble();
            *M_blockPreconditioner += *matrixNoBC;
Alessio Fumagalli's avatar
Alessio Fumagalli committed
1611
1612
            matrix_Type tempMatrix ( *matrixNoBC );
            matrixNoBC.reset ( new matrix_Type ( M_localMap, tempMatrix.meanNumEntries() ) );