NavierStokesSolverPC.hpp 20.2 KB
Newer Older
prudhomm's avatar
prudhomm committed
1
2
/*
  This file is part of the LifeV library
3
  Copyright (C) 2001,2002,2003,2004 EPFL, INRIA and Politechnico di Milano
prudhomm's avatar
prudhomm committed
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
  
  This library 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 2.1 of the License, or (at your option) any later version.
  
  This library 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, write to the Free Software
  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
*/
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
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
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
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
/*!
  \file NavierStokesSolverPC.h
  \author M.A. Fernandez and A. Gauthier
  \date 11/2002 
  \version 1.0
  \author A. Veneziani
  \date 05/2003
  \version 1.1


  \brief This file contains a NavierStokes solver class which implements a semi-implicit scheme with an exact factorization. 
         Preconditioning of the Schur Complement is done by an algebraic Chorin-Temam pressure-corrected preconditioner 
         Added a class for handling high order discretization schemes (A. Veneziani) 
*/

#ifndef _NAVIERSTOKESSOLVERPC_H_
#define _NAVIERSTOKESSOLVERPC_H_

#include "NavierStokesHandler.hpp"
#include "elemMat.hpp"
#include "elemVec.hpp"
#include "elemOper.hpp"
#include "values.hpp"
#include "pattern.hpp"
#include "assemb.hpp"
#include "bc_manage.hpp"
#include "algebraic_facto.hpp"
#include "bcCond.hpp"
#include "chrono.hpp"
#include "dataAztec.hpp"
#include "bdfNS.hpp"
#include "openDX_wrtrs.hpp"
#include <string>

/*! 
  \class NavierStokesSolverPC

   This class implements an NavierStokes solver via exact factorization. Preconditioning of the 
  Schur Complement is done by an algebraic Chorin-Temam pressure-corrected preconditioner 

*/
template<typename Mesh>
class NavierStokesSolverPC:
public NavierStokesHandler<Mesh> {
 
 public:
  
  typedef  typename  NavierStokesHandler<Mesh>::Function Function; 

  //! Constructor 
  /*!
    \param data_file GetPot data file
    \param refFE_u reference FE for the velocity
    \param refFE_p reference FE for the pressure
    \param Qr_u volumic quadrature rule for the velocity
    \param bdQr_u surface quadrature rule for the velocity 
    \param Qr_p volumic quadrature rule for the pressure
    \param bdQr_p surface quadrature rule for the pressure
    \param BCh_u boundary conditions for the velocity
    \param ord_bdf order of the Bdf time advancing scheme + incremental for the pressure
  */
  NavierStokesSolverPC(const GetPot& data_file, const RefFE& refFE_u, const RefFE& refFE_p, const QuadRule& Qr_u,
	    const QuadRule& bdQr_u, const QuadRule& Qr_p, const QuadRule& bdQr_p, BC_Handler& BCh_u);

  //! Update the right  hand side  for time advancing 
  /*! 
    \param source volumic source  
    \param time present time
  */
  void timeAdvance(const Function source, const Real& time);

  //! Update convective term, bc treatment and solve the linearized ns system
  void iterate(const Real& time);

  // ! Residual Computation
  PhysVectUnknown<Vector> residual();

  // ! Shear stress computation ***** Prova Agosto 2003
  void ShearStressCompute(string filename_sstress, string fe_type);


 private:

  //! Block pattern of C: rho/dt*Vmass + mu*Vstiff operator
  MSRPatt _pattC_block;
  
  //! Pattern for C
  MixedPattern<3,3,MSRPatt> _pattC;

  //! Block pattern of D: Vdiv operator
  CSRPatt _pattD_block;

  //! Pattern for D
  MixedPattern<1,3,CSRPatt> _pattD;

  //! Block  pattern of trD: transpose Vdiv operator trVdiv
  CSRPatt _pattDtr_block;
  
  //! Pattern  of trD
  MixedPattern<3,1,CSRPatt> _pattDtr;

  //! Matrix D: Vdiv operator
  MixedMatr<1,3,CSRPatt,double> _D;
 
  //! Matrix trD: transpose Vdiv operator trVdiv
  MixedMatr<3,1,CSRPatt,double> _trD;

  //! Matrix HinvDtr:  H^{-1}D^T
  MixedMatr<3,1,CSRPatt,double> _HinvDtr;

  //! Matrix M_u: Vmass
  MixedMatr<3,3,MSRPatt,double> _M_u;
 
  //! Matrix HinvC: H^{-1}C
  MixedMatr<3,3,MSRPatt,double> _HinvC;
 
  //! Matrix C: rho/dt*Vmass + mu*Vstiff operator
  MixedMatr<3,3,MSRPatt,double> _CStokes;
  
  //! Matrix C: rho/dt*Vmass + mu*Vstiff operator + Convective_term
  MixedMatr<3,3,MSRPatt,double> _C;

 
  //! H diag matrix: H= diag( _M_u )/sum( diag( _M_u ) ) where _M_u = mass * rho / dt
  vector<double>  _H;

  //! Elementary matrices and vectors
  ElemMat _elmatC; //velocity stiffnes 
  ElemMat _elmatM_u; //velocity mass
  ElemMat _elmatDtr; // vel_i * pres_j
  ElemVec _elvec; // Elementary right hand side

  //! Right  hand  side for the velocity
  PhysVectUnknown<Vector> _f_u;
  
  //!  This vector contains the product C^{-1}*trD*P where P is the pressure, solution
  //!  of the system (ii).
  PhysVectUnknown<Vector> _invCtrDP;

  // A Veneziani August 2003 ****************
  //! Matrix Cnobc: rho/dt*Vmass + mu*Vstiff operator + Convective_term WITHOUT BC (for the residual computation)
  MixedMatr<3,3,MSRPatt,double> _CnoBc;

  //! Matrix trD: transpose Vdiv operator trVdiv WITHOUT BC (for the residual computation)
  MixedMatr<3,1,CSRPatt,double> _trDnoBc;

  //! Right  hand  side for the velocity WITHOUT BC
  PhysVectUnknown<Vector> _f_u_noBc;
  // REM This solution is just to start: Miguel suggested a different way,
  // dealing only with the submatrices involved in the Dirichlet elimination.
  // This is a good idea....still to be done


  DataAztec _dataAztec_i;
  DataAztec _dataAztec_ii;
  DataAztec _dataAztec_s;

  //! DataFactorisation: data passed to matrix-vector product are stored in the class
  DataFactorisation<
    MixedMatr<3,3,MSRPatt,double>,
    MixedMatr<1,3,CSRPatt,double>,
    MixedMatr<3,1,CSRPatt,double>,
    vector<double>,
    MSRMatr<double>,
    Vector> _factor_data;
};


//
//                                         IMPLEMENTATION
//
template<typename Mesh> NavierStokesSolverPC<Mesh>::
NavierStokesSolverPC(const GetPot& data_file, const RefFE& refFE_u, const RefFE& refFE_p, const QuadRule& Qr_u,
	  const QuadRule& bdQr_u, const QuadRule& Qr_p, const QuadRule& bdQr_p, BC_Handler& BCh_u):
  NavierStokesHandler<Mesh>(data_file,refFE_u,refFE_p,Qr_u,bdQr_u,Qr_p,bdQr_p,BCh_u),
     _pattC_block(_dof_u),
     _pattC(_pattC_block,"diag"),
     _pattD_block(_dof_p,_dof_u),
     _pattD(_pattD_block),
     _pattDtr_block(_dof_u,_dof_p),
     _pattDtr(_pattDtr_block),
     _D(_pattD),
     _trD(_pattDtr),
     _HinvDtr(_pattDtr),
     _M_u(_pattC),
     _HinvC(_pattC),
     _CStokes(_pattC),
     _C(_pattC),
     _H(_pattC.nRows()), 
     _elmatC(_fe_u.nbNode,nDimensions,nDimensions), 
     _elmatM_u(_fe_u.nbNode,nDimensions,nDimensions),
     _elmatDtr(_fe_u.nbNode,nDimensions,0,_fe_p.nbNode,0,1), 
     _elvec(_fe_u.nbNode,nDimensions), 
     _f_u(_dim_u),
     _invCtrDP(_dim_u),
      _CnoBc(_pattC),_trDnoBc(_pattDtr),_f_u_noBc(_dim_u),
     _dataAztec_i(data_file,"fluid/aztec_i"),
     _dataAztec_ii(data_file,"fluid/aztec_ii"),
     _dataAztec_s(data_file,"fluid/aztec_s"),
     _factor_data(_C,_D,_trD,_H,_HinvC,_HinvDtr,_invCtrDP.vec(),_dataAztec_i,_dataAztec_s,_BCh_u.fullEssential()) {
  
  cout << endl;
  cout << "O-  Pressure unknowns: " << _dim_p     << endl; 
  cout << "O-  Velocity unknowns: " << _dim_u     << endl<<endl;
  cout << "O-  Computing mass and Stokes matrices... ";  
    
  Chrono chrono;
  chrono.start();

  // Matrices initialization 
  _D.zeros();
  _trD.zeros();
  _M_u.zeros();
  _CStokes.zeros();
  _C.zeros();
  
  _CnoBc.zeros();
  _trDnoBc.zeros();
  

  // Number of velocity components  
  UInt nc_u=_u.nbcomp();
  
  //inverse of dt:
  Real dti=1./_dt;


  // *******************************************************
  // Coefficient of the mass term at time t^{n+1}
  Real first_coeff = _bdf.bdf_u().coeff_der(0);
  cout << endl;
  cout << "Bdf NS first coeff " << first_coeff << endl; 

  _bdf.bdf_u().showMe();
  _bdf.bdf_p().showMe();

  // Auxiliary matrix
  ElemMat elmatM_u_St(_fe_u.nbNode,nDimensions,nDimensions);

  // Elementary computation and matrix assembling  
  // Loop on elements
  for(UInt i = 1; i <= _mesh.numVolumes(); i++){

    _fe_p.update( _mesh.volumeList(i) ); // just to provide the id number in the assem_mat_mixed
    _fe_u.updateFirstDerivQuadPt(_mesh.volumeList(i));
    
    _elmatC.zero();
    _elmatM_u.zero();
     elmatM_u_St.zero();
    _elmatDtr.zero();

    stiff(_mu,_elmatC,_fe_u,0,0,nDimensions); 
  // *******************************************************
    mass(first_coeff*_rho*dti,elmatM_u_St,_fe_u,0,0,nDimensions);
    mass(_rho*dti,_elmatM_u,_fe_u,0,0,nDimensions);
    // stiffness + mass
    
    _elmatC.mat() += elmatM_u_St.mat();
    
    for(UInt ic=0;ic<nc_u;ic++){
    
      // stiffness
      assemb_mat(_CStokes,_elmatC,_fe_u,_dof_u,ic,ic);
      
      // mass
      assemb_mat(_M_u,_elmatM_u,_fe_u,_dof_u,ic,ic);
     
      // div
      grad(ic,1.0,_elmatDtr,_fe_u,_fe_p,ic,0);
     
      // assembling
      assemb_mat_mixed(_trD,_elmatDtr,_fe_u,_fe_p,_dof_u,_dof_p,ic,0);
      assemb_tr_mat_mixed(1.0,_D,_elmatDtr,_fe_p,_fe_u,_dof_p,_dof_u,0,ic); 
    }
  }

  // H diag matrix: H= diag( _M_u )/sum( diag( _M_u ) ) where _M_u = mass * rho / dt * first_coeff_bdf_scheme
  _H=_M_u.giveDiag();
  double sum=accumulate(_H.begin(),_H.end(),0.0);
  // *******************************************************
  // Matrix equilibration: 
  for (UInt i=0; i<_H.size();i++) 
    _H[i]=_H[i]*first_coeff*dti/sum; // H = lumping of first_coeff/dt*Mass
    
  chrono.stop();
  cout << "done in " << chrono.diff() << " s." << endl;
}

template<typename Mesh>  
void NavierStokesSolverPC<Mesh>::
timeAdvance(const Function source, const Real& time) {

  cout << endl;
  cout << "O== Now we are at time "<< time << " s." << endl;

  // Number of velocity components  
  UInt nc_u=_u.nbcomp();

  cout << "  o-  Updating mass term on right hand side... ";


  Chrono chrono;
  chrono.start();

  // Right hand side for the velocity at time
  _f_u.vec()=0.;

  // loop on volumes: assembling source term
  for(UInt i=1; i<=_mesh.numVolumes(); ++i){
    _elvec.zero();
    _fe_u.update(_mesh.volumeList(i));

    for (UInt ic=0; ic<nc_u; ++ic){ 
      compute_vec(source,_elvec,_fe_u,time,ic); // compute local vector
      assemb_vec(_f_u.vec(),_elvec,_fe_u,_dof_u,ic); // assemble local vector into global one       
    }
  }

  // ******************************************************* 
  _f_u.vec() += _M_u*_bdf.bdf_u().time_der(); //_M_u is the mass matrix divided by the time step
  //  _f_u.vec() += _M_u * _u.vec();
  chrono.stop();
  cout << "done in " << chrono.diff() << " s." << endl;
}


template<typename Mesh>  
void NavierStokesSolverPC<Mesh>::
iterate(const Real& time) {

  // Number of velocity components  
  UInt nc_u=_u.nbcomp();
  Vector u_extrap = _bdf.bdf_u().extrap();

  Chrono  chrono;

  // C = CStokes + convective term
  chrono.start();
  _C=_CStokes;
  chrono.stop();


  cout << "  o-  Stokes matrix was copied in " << chrono.diff() << "s." << endl;
  cout << "  o-  Updating convective term... ";
 
  chrono.start();

  // loop on volumes
  for(UInt i=1; i<=_mesh.numVolumes(); ++i){
  
    _fe_u.updateFirstDeriv(_mesh.volumeList(i)); // as updateFirstDer

    _elmatC.zero();

    UInt eleID = _fe_u.currentId();
    // Non linear term, Semi-implicit approach      
    // ULoc contains the velocity values in the nodes
    for (UInt k=0 ; k<(UInt)_fe_u.nbNode ; k++){
	UInt  iloc = _fe_u.patternFirst(k);
	for (UInt ic=0; ic<nc_u; ++ic){     
	  UInt ig=_dof_u.localToGlobal(eleID,iloc+1)-1+ic*_dim_u;       
	  _elvec.vec()[iloc+ic*_fe_u.nbNode] = u_extrap(ig); 
	}
    }

    // loop on components 
    for (UInt ic=0; ic<nc_u; ++ic){
      // compute local convective term and assembling     	
      grad(0,_elvec,_elmatC,_fe_u,_fe_u,ic,ic);
      grad(1,_elvec,_elmatC,_fe_u,_fe_u,ic,ic);
      grad(2,_elvec,_elmatC,_fe_u,_fe_u,ic,ic);
      assemb_mat(_C,_elmatC,_fe_u,_dof_u,ic,ic);
    }
  }
  chrono.stop();
  cout << "done in " << chrono.diff() << "s." << endl;

  // QUI VENGONO APPLICATE LE BC
  _CnoBc=_C;
  _trDnoBc=_trD;
  _f_u_noBc.vec()=_f_u.vec(); 
  //for (UInt myindex=0;myindex<_dim_u;myindex++) _f_u_noBc.vec()[myindex] = _f_u.vec()[myindex];
 
  _f_u.vec() -= _trD*_bdf.bdf_p().extrap(); // INCREMENTAL correction for the pressure: IT MUST BE AFTER THE INITIALIZATION of _f_u_noBC

  // for BC treatment (done at each time-step)
  Real tgv=1.e02; 

  cout << "  o-  Applying boundary conditions... ";
  chrono.start(); 
  // BC manage for the velocity
  if ( !_BCh_u.bdUpdateDone() )  
    _BCh_u.bdUpdate(_mesh, _feBd_u, _dof_u);
  bc_manage(_C, _trD, _f_u.vec(), _mesh, _dof_u, _BCh_u, _feBd_u, tgv, time);
  chrono.stop();
  cout << "done in " << chrono.diff() << "s." << endl;
 
  //matrices HinvDtr:
  MultInvDiag(_H, _trD, _HinvDtr);
 
  // AZTEC specifications for each system
  int    data_org_i1[AZ_COMM_SIZE];   // data organisation for C1
  int    data_org_i2[AZ_COMM_SIZE];   // data organisation for C2
  int    data_org_i3[AZ_COMM_SIZE];   // data organisation for C3
 
  // identical for each system
  int    proc_config_i[AZ_PROC_SIZE];// Processor information:  
  int    options_i[AZ_OPTIONS_SIZE]; // Array used to select solver options.
  double params_i[AZ_PARAMS_SIZE];   // User selected solver paramters.
  double status_i[AZ_STATUS_SIZE];   // Information returned from AZ_solve()
                                     // indicating success or failure.
  AZ_set_proc_config(proc_config_i, AZ_NOT_MPI);

  //AZTEC matrix and preconditioner
  AZ_MATRIX *C1, *C2, *C3;
  AZ_PRECOND *prec_C1, *prec_C2, *prec_C3;

  int N_eq_i= _dim_u; // number of DOF for each component

  // set each block
  C1= AZ_matrix_create(N_eq_i);
  C2= AZ_matrix_create(N_eq_i);
  C3= AZ_matrix_create(N_eq_i);

  // create preconditioner for each block
  prec_C1= AZ_precond_create(C1, AZ_precondition, NULL);
  prec_C2= AZ_precond_create(C2, AZ_precondition, NULL);
  prec_C3= AZ_precond_create(C3, AZ_precondition, NULL);

  // data_org assigned "by hands" while no parallel computation is performed
  data_org_i1[AZ_N_internal]= N_eq_i;
  data_org_i1[AZ_N_border]= 0;
  data_org_i1[AZ_N_external]= 0;
  data_org_i1[AZ_N_neigh]= 0;
  data_org_i1[AZ_name]= DATA_NAME_AZTEC1;
  
  data_org_i2[AZ_N_internal]= N_eq_i;
  data_org_i2[AZ_N_border]= 0;
  data_org_i2[AZ_N_external]= 0;
  data_org_i2[AZ_N_neigh]= 0;
  data_org_i2[AZ_name]= DATA_NAME_AZTEC2;

  data_org_i3[AZ_N_internal]= N_eq_i;
  data_org_i3[AZ_N_border]= 0;
  data_org_i3[AZ_N_external]= 0;
  data_org_i3[AZ_N_neigh]= 0;
  data_org_i3[AZ_name]= DATA_NAME_AZTEC3;

  // --------------------------------------
  // (ii) (D*C^{-1}*trD) * \delta P =  D*V
  // --------------------------------------
  
  // AZTEC specifications for the second system
  int    proc_config_ii[AZ_PROC_SIZE];// Processor information:
  //  proc_config[AZ_node] = node name
  //  proc_config[AZ_N_procs] = # of nodes
  int    options_ii[AZ_OPTIONS_SIZE]; // Array used to select solver options.
  double params_ii[AZ_PARAMS_SIZE];   // User selected solver paramters.
  double status_ii[AZ_STATUS_SIZE];   // Information returned from AZ_solve()
                                     // indicating success or failure.
  //

  AZ_set_proc_config(proc_config_ii, AZ_NOT_MPI);

  //AZTEC matrix for A_ii=(D*C^{-1}*trD)
  AZ_MATRIX *A_ii;
  AZ_PRECOND *pILU_ii;

  int N_eq_ii= _p.size();

  A_ii= AZ_matrix_create(N_eq_ii);
  // data containing the matrices C, D, trD and H as pointers
  // are passed through A_ii and pILU_ii:
  AZ_set_MATFREE(A_ii, &_factor_data,
		 my_matvec_block<
		 MixedMatr<3,3,MSRPatt,double>,
		 MixedMatr<1,3,CSRPatt,double>,
		 MixedMatr<3,1,CSRPatt,double>,
		 vector<double>,
		 MSRMatr<double>,
		 Vector>);
  
  pILU_ii= AZ_precond_create(A_ii, my_precSchur_PC<
			     MixedMatr<3,3,MSRPatt,double>,
			     MixedMatr<1,3,CSRPatt,double>,
			     MixedMatr<3,1,CSRPatt,double>,
			     vector<double>,
			     MSRMatr<double>,
			     Vector>, &_factor_data);
  

  _dataAztec_ii.aztecOptionsFromDataFile(options_ii,params_ii);
   
  // user preconditioning:
  options_ii[AZ_precond]  = AZ_user_precond;
  
  // RHS of the linear system (ii)
  Vector vec_DV(_p.size());

  //matrices HinvC (depends on time):
  MultInvDiag(_H, _C, _HinvC);

  // ---------------
  // (i) C * V = F_V
  // ---------------
  cout << "  o-  Solving first system... "; 

  AZ_set_MSR(C1, (int*) _pattC_block.giveRaw_bindx(),
	     (double*) _C.giveRaw_value(0,0),
	     data_org_i1, 0, NULL, AZ_LOCAL);
  AZ_set_MSR(C2, (int*) _pattC_block.giveRaw_bindx(),
	     (double*) _C.giveRaw_value(1,1),
	     data_org_i2, 0, NULL, AZ_LOCAL);
  AZ_set_MSR(C3, (int*) _pattC_block.giveRaw_bindx(),
	     (double*) _C.giveRaw_value(2,2),
	     data_org_i3, 0, NULL, AZ_LOCAL);

  _dataAztec_i.aztecOptionsFromDataFile(options_i,params_i);
  
  //keep C factorisation and preconditioner reused in my_matvec
  options_i[AZ_keep_info]= 1;               // keep information
   
  // ---------------
  // (i) C * V = F_V (= forcing term + bc - D^T*P^n)
  // ---------------
  // intermediate velocity computation
  
  chrono.start();

  //for each block
  AZ_iterate(_u.giveVec(), _f_u.giveVec(), options_i,params_i, status_i,
	     proc_config_i, C1, prec_C1, NULL);
  AZ_iterate(_u.giveVec()+_dim_u, _f_u.giveVec()+_dim_u,  options_i,params_i,
	     status_i, proc_config_i, C2, prec_C2, NULL);
  AZ_iterate(_u.giveVec()+2*_dim_u, _f_u.giveVec()+2*_dim_u,  options_i,params_i,
	     status_i, proc_config_i, C3, prec_C3, NULL);

  //
  chrono.stop();
  cout << "done in " << chrono.diff() << " s." << endl;

  // ---------------------------------------------------
  // (ii) (D*C^(-1)*trD) * \delta P = D*C^{-1}*F_V = D*V
  // ---------------------------------------------------

  // RHS of the linear system (ii)
  vec_DV = _D*_u.vec();
  _p.vec() = 0.0; // AT this point, this vector stands for the "pressure increment"
  
  // case of pure Dirichlet BCs:
  if (_BCh_u.fullEssential()) {
    vec_DV[_dim_p-1]   = 0.0; // correction of the right hand side.
    _p.vec()[_dim_p-1] = 0.0; // pressure value at the last node.
  }

  cout << "  o-  Solving second system... ";

  chrono.start();
  AZ_iterate(_p.giveVec(), &vec_DV[0],  options_ii,params_ii, status_ii,
	     proc_config_ii, A_ii, pILU_ii, NULL);

  chrono.stop();
  cout << "done in " << chrono.diff() << " s." << endl;
  
  // ------------------------------------
  // (iii) V = V-(C^(-1)*trD) * \delta P
  // ------------------------------------

  
  // everything is done...
  _u.vec() = _u.vec() - _invCtrDP.vec();
  cout << "  o-  Velocity updated" << endl;

  // ******************************************************* 
  // This is the REAL pressure (not the increment)
  _p.vec() += _bdf.bdf_p().extrap();
  cout << "  o-  Pressure updated" << endl;

  // *******************************************************
  // update the array of the previous solutions
  _bdf.bdf_u().shift_right(_u.vec());
  _bdf.bdf_p().shift_right(_p.vec());
  

  // destroy Ci and prec_Ci
  AZ_matrix_destroy(&C1);
  AZ_precond_destroy(&prec_C1);
  AZ_matrix_destroy(&C2);
  AZ_precond_destroy(&prec_C2);
  AZ_matrix_destroy(&C3);
  AZ_precond_destroy(&prec_C3);
  // destroy A_ii and pILU_ii
  AZ_matrix_destroy(&A_ii);
  AZ_precond_destroy(&pILU_ii);
}


////////////////

template<typename Mesh>  
PhysVectUnknown<Vector> NavierStokesSolverPC<Mesh>::residual()
 {
  Chrono chrono;
  PhysVectUnknown<Vector> r(_dim_u);
  cout << "  o- Computing the residual...";
  chrono.start();
  r.vec() = _f_u_noBc.vec()-_CnoBc*_u.vec() - _trDnoBc*_p.vec();
  chrono.stop();
  cout << "done in " << chrono.diff() << "s" << endl;
  return r;
  }

template<typename Mesh>  
void NavierStokesSolverPC<Mesh>::ShearStressCompute(string filename_stress, string fe_type)
{
  Vector residual(this->residual().vec());
  UInt ss = residual.size()/NDIM;
  Vector sstress(_ns_post_proc.compute_sstress(residual,(UInt)NDIM));
  // just a stupid way for writing the shear stress in OpenDx or Medit formats,
  // exploting the existent subroutines
  UInt s=sstress.size()/NDIM;
  UInt where;
  for (UInt i=0;i<s;i++){
    where=_ns_post_proc.fBdToIn()[i];       
    //    cout << residual.size() << " " << s << " " << where << " " << i << endl;
   for (UInt j=0;j<NDIM;j++) residual[where-1+j*ss]=sstress[i+j*s];
  }  
  wr_opendx_header(filename_stress, _mesh, _dof_u,_fe_u,fe_type);
  wr_opendx_vector(filename_stress,"ShearStress",residual,NDIM);
}

 

#endif