File: predvv.f90

package info (click to toggle)
code-saturne 4.3.3%2Brepack-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 77,992 kB
  • sloc: ansic: 281,257; f90: 122,305; python: 56,490; makefile: 3,915; xml: 3,285; cpp: 3,183; sh: 1,139; lex: 176; yacc: 101; sed: 16
file content (1947 lines) | stat: -rw-r--r-- 61,294 bytes parent folder | download
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
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
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
713
714
715
716
717
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
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
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
824
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
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
954
955
956
957
958
959
960
961
962
963
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
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
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
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
!-------------------------------------------------------------------------------

! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2016 EDF S.A.
!
! This program is free software; you can redistribute it and/or modify it under
! the terms of the GNU General Public License as published by the Free Software
! Foundation; either version 2 of the License, or (at your option) any later
! version.
!
! This program 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 General Public License for more
! details.
!
! You should have received a copy of the GNU General Public License along with
! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
! Street, Fifth Floor, Boston, MA 02110-1301, USA.

!-------------------------------------------------------------------------------

!===============================================================================
! Function:
! ---------

!> \file predvv.f90
!>
!> \brief This subroutine performs the velocity prediction step of the Navier
!> Stokes equations for incompressible or slightly compressible flows for
!> the coupled velocity components solver.
!>
!> - at the first call, the predicted velocities are computed and also
!>   an estimator on the predicted velocity is computed.
!>
!> - at the second call, a global estimator on Navier Stokes is computed.
!>   This second call is done after the correction step (\ref resopv).
!>
!-------------------------------------------------------------------------------

!-------------------------------------------------------------------------------
! Arguments
!______________________________________________________________________________.
!  mode           name          role                                           !
!______________________________________________________________________________!
!> \param[in]     iappel        call number (1 or 2)
!> \param[in]     nvar          total number of variables
!> \param[in]     nscal         total number of scalars
!> \param[in]     iterns        index of the iteration on Navier-Stokes
!> \param[in]     ncepdp        number of cells with head loss
!> \param[in]     ncesmp        number of cells with mass source term
!> \param[in]     nfbpcd        number of faces with condensation source terms
!> \param[in]     ncmast        number of cells with condensation source terms
!> \param[in]     icepdc        index of cells with head loss
!> \param[in]     icetsm        index of cells with mass source term
!> \param[in]     ifbpcd        index of faces with condensation source terms
!> \param[in]     ltmast        index of cells with condensation source terms
!> \param[in]     itypsm        type of mass source term for the variables
!> \param[in]     dt            time step (per cell)
!> \param[in]     vel           velocity
!> \param[in]     vela          velocity at the previous time step
!> \param[in]     flumas        internal mass flux (depending on iappel)
!> \param[in]     flumab        boundary mass flux (depending on iappel)
!> \param[in]     tslagr        coupling term for the Lagrangian module
!> \param[in]     coefav        boundary condition array for the variable
!>                               (explicit part)
!> \param[in]     coefbv        boundary condition array for the variable
!>                               (implicit part)
!> \param[in]     cofafv        boundary condition array for the diffusion
!>                               of the variable (explicit part)
!> \param[in]     cofbfv        boundary condition array for the diffusion
!>                               of the variable (implicit part)
!> \param[in]     ckupdc        work array for the head loss
!> \param[in]     smacel        variable value associated to the mass source
!>                               term (for ivar=ipr, smacel is the mass flux
!>                               \f$ \Gamma^n \f$)
!> \param[in]     spcond        variable value associated to the condensation
!>                              source term (for ivar=ipr, spcond is the flow rate
!>                              \f$ \Gamma_{s, cond}^n \f$)
!> \param[in]     svcond        variable value associated to the condensation
!>                              source term (for ivar=ipr, svcond is the flow rate
!>                              \f$ \Gamma_{v, cond}^n \f$)
!> \param[in]     frcxt         external forces making hydrostatic pressure
!> \param[in]     trava         working array for the velocity-pressure coupling
!> \param[in]     ximpa         same
!> \param[in]     uvwk          same (stores the velocity at the previous iteration)
!> \param[in]     dfrcxt        variation of the external forces
!                               making the hydrostatic pressure
!> \param[in]     grdphd        hydrostatic pressure gradient to handle the
!>                              imbalance between the pressure gradient and
!>                              gravity source term
!> \param[in]     tpucou        non scalar time step in case of
!>                              velocity pressure coupling
!> \param[in]     trav          right hand side for the normalizing
!>                              the residual
!> \param[in]     viscf         visc*surface/dist aux faces internes
!> \param[in]     viscb         visc*surface/dist aux faces de bord
!> \param[in]     viscfi        same as viscf for increments
!> \param[in]     viscbi        same as viscb for increments
!> \param[in]     secvif        secondary viscosity at interior faces
!> \param[in]     secvib        secondary viscosity at boundary faces
!> \param[in]     w1            working array
!> \param[in]     w7            working array
!> \param[in]     w8            working array
!> \param[in]     w9            working array
!_______________________________________________________________________________

subroutine predvv &
 ( iappel ,                                                       &
   nvar   , nscal  , iterns ,                                     &
   ncepdp , ncesmp , nfbpcd , ncmast ,                            &
   icepdc , icetsm , ifbpcd , ltmast ,                            &
   itypsm ,                                                       &
   dt     , vel    , vela   ,                                     &
   flumas , flumab ,                                              &
   tslagr , coefav , coefbv , cofafv , cofbfv ,                   &
   ckupdc , smacel , spcond , svcond , frcxt  , grdphd ,          &
   trava  , ximpa  , uvwk   , dfrcxt , tpucou , trav   ,          &
   viscf  , viscb  , viscfi , viscbi , secvif , secvib ,          &
   w1     , w7     , w8     , w9     )

!===============================================================================

!===============================================================================
! Module files
!===============================================================================

use paramx
use dimens, only: ndimfb
use numvar
use entsor
use cstphy
use cstnum
use optcal
use parall
use period
use lagran
use ppppar
use ppthch
use ppincl
use cplsat
use ihmpre, only: iihmpr
use mesh
use rotation
use turbomachinery
use cs_f_interfaces
use cs_c_bindings
use cfpoin
use field
use field_operator
use pointe, only: gamcav
use cavitation
use cs_tagms, only:s_metal

!===============================================================================

implicit none

! Arguments

integer          iappel
integer          nvar   , nscal  , iterns
integer          ncepdp , ncesmp , nfbpcd , ncmast

integer          icepdc(ncepdp)
integer          icetsm(ncesmp), itypsm(ncesmp,nvar)
integer          ifbpcd(nfbpcd)
integer          ltmast(ncelet)

double precision dt(ncelet)
double precision flumas(nfac), flumab(nfabor)
double precision tslagr(ncelet,*)
double precision ckupdc(ncepdp,6), smacel(ncesmp,nvar)
double precision spcond(nfbpcd,nvar), svcond(ncelet,nvar)
double precision frcxt(3,ncelet), dfrcxt(3,ncelet)
double precision grdphd(3, ncelet)
double precision trava(ndim,ncelet)
double precision ximpa(ndim,ndim,ncelet),uvwk(ndim,ncelet)
double precision tpucou(6, ncelet)
double precision trav(3,ncelet)
double precision viscf(*), viscb(nfabor)
double precision viscfi(*), viscbi(nfabor)
double precision secvif(nfac), secvib(nfabor)
double precision w1(ncelet)
double precision w7(ncelet), w8(ncelet), w9(ncelet)
double precision coefav(3  ,ndimfb)
double precision cofafv(3  ,ndimfb)
double precision coefbv(3,3,ndimfb)
double precision cofbfv(3,3,ndimfb)

double precision vel   (3  ,ncelet)
double precision vela  (3  ,ncelet)

! Local variables

integer          f_id  , iel   , ielpdc, ifac  , isou  , itypfl, n_fans
integer          iccocg, inc   , iprev , init  , ii    , jj
integer          nswrgp, imligp, iwarnp
integer          iswdyp, idftnp
integer          iconvp, idiffp, ndircp, nswrsp
integer          ircflp, ischcp, isstpp, iescap
integer          iflmb0, nswrp
integer          idtva0, icvflb
integer          jsou  , ivisep, imasac
integer          ivoid(1)

double precision rnorm , vitnor
double precision romvom, drom  , rom
double precision epsrgp, climgp, extrap, relaxp, blencp, epsilp
double precision epsrsp
double precision vit1  , vit2  , vit3, xkb, pip, pfac, pfac1
double precision cpdc11, cpdc22, cpdc33, cpdc12, cpdc13, cpdc23
double precision d2s3  , thetap, thetp1, thets , dtsrom
double precision diipbx, diipby, diipbz
double precision ccorio
double precision dvol

double precision rvoid(1)

! Working arrays
double precision, allocatable, dimension(:,:) :: eswork
double precision, allocatable, dimension(:,:) :: grad
double precision, dimension(:,:), allocatable :: smbr
double precision, dimension(:,:,:), allocatable :: fimp
double precision, dimension(:,:), allocatable :: gavinj
double precision, dimension(:,:), allocatable :: tsexp
double precision, dimension(:,:,:), allocatable :: tsimp
double precision, allocatable, dimension(:,:) :: viscce
double precision, dimension(:,:), allocatable :: vect
double precision, dimension(:), allocatable :: xinvro
double precision, dimension(:), pointer :: brom, crom, croma, pcrom
double precision, dimension(:), pointer :: coefa_k, coefb_k
double precision, dimension(:), pointer :: coefa_p, coefb_p
double precision, dimension(:,:), allocatable :: rij
double precision, dimension(:), pointer :: coef1, coef2, coef3, coef4, coef5, coef6
double precision, dimension(:,:), pointer :: coefap
double precision, dimension(:,:,:), pointer :: coefbp
double precision, dimension(:,:), allocatable :: coefat
double precision, dimension(:,:,:), allocatable :: coefbt
double precision, dimension(:,:), allocatable :: tflmas, tflmab
double precision, dimension(:,:), allocatable :: divt
double precision, dimension(:), allocatable ::xnormp
double precision, dimension(:,:), pointer :: forbr, c_st_vel
double precision, dimension(:), pointer :: cvara_pr, cvara_k
double precision, dimension(:), pointer :: cvara_r11, cvara_r22, cvara_r33
double precision, dimension(:), pointer :: cvara_r12, cvara_r23, cvara_r13
double precision, dimension(:,:), pointer :: cvara_rij
double precision, dimension(:), pointer :: viscl, visct, c_estim
double precision, allocatable, dimension(:) :: surfbm
double precision, dimension(:,:), pointer :: lapla
double precision, dimension(:), pointer :: cpro_tsrho

!===============================================================================

!===============================================================================
! 1. Initialization
!===============================================================================

! Allocate temporary arrays
allocate(smbr(3,ncelet))
allocate(fimp(3,3,ncelet))
allocate(tsexp(3,ncelet))
allocate(tsimp(3,3,ncelet))
if (idften(iu).eq.6) allocate(viscce(6,ncelet))

! Allocate a temporary array for the prediction-stage error estimator
if (iescal(iespre).gt.0) then
  allocate(eswork(3,ncelet))
endif

! Reperage de rho au bord
call field_get_val_s(ibrom, brom)
! Reperage de rho courant (ie en cas d'extrapolation rho^n+1/2)
call field_get_val_s(icrom, crom)
! Reperage de rho^n en cas d'extrapolation
if (iroext.gt.0.or.idilat.gt.1) then
  call field_get_val_prev_s(icrom, croma)
endif


if (iappel.eq.2) then
  if (iforbr.ge.0 .and. iterns.eq.1 .or. icavit.ge.0) then
    call field_get_val_s(ivarfl(ipr), cvara_pr)
  endif
  if(iforbr.ge.0 .and. iterns.eq.1                                          &
     .and. (itytur.eq.2 .or. itytur.eq.5 .or. iturb.eq.60) .and. igrhok.eq.1) then
    call field_get_val_s(ivarfl(ik), cvara_k)
  endif
  if (itytur.eq.3.and.iterns.eq.1) then
    if (irijco.eq.1) then
      call field_get_val_v(ivarfl(irij), cvara_rij)
    else
      call field_get_val_s(ivarfl(ir11), cvara_r11)
      call field_get_val_s(ivarfl(ir22), cvara_r22)
      call field_get_val_s(ivarfl(ir33), cvara_r33)
      call field_get_val_s(ivarfl(ir12), cvara_r12)
      call field_get_val_s(ivarfl(ir23), cvara_r23)
      call field_get_val_s(ivarfl(ir13), cvara_r13)
    endif
  endif
else
  if (iforbr.ge.0 .and. iterns.eq.1 .or. icavit.ge.0) then
    call field_get_val_prev_s(ivarfl(ipr), cvara_pr)
  endif
  if(iforbr.ge.0 .and. iterns.eq.1                                          &
     .and. (itytur.eq.2 .or. itytur.eq.5 .or. iturb.eq.60) .and. igrhok.eq.1) then
      call field_get_val_prev_s(ivarfl(ik), cvara_k)
  endif
  if (itytur.eq.3.and.iterns.eq.1) then
    if (irijco.eq.1) then
      call field_get_val_v(ivarfl(irij), cvara_rij)
    else
      call field_get_val_s(ivarfl(ir11), cvara_r11)
      call field_get_val_s(ivarfl(ir22), cvara_r22)
      call field_get_val_s(ivarfl(ir33), cvara_r33)
      call field_get_val_s(ivarfl(ir12), cvara_r12)
      call field_get_val_s(ivarfl(ir23), cvara_r23)
      call field_get_val_s(ivarfl(ir13), cvara_r13)
    endif
  endif
endif

if (iforbr.ge.0 .and. iterns.eq.1) call field_get_val_v(iforbr, forbr)

! Theta relatif aux termes sources explicites
thets  = thetsn
if (isno2t.gt.0) then
  call field_get_key_int(ivarfl(iu), kstprv, f_id)
  call field_get_val_v(f_id, c_st_vel)
else
  c_st_vel => null()
endif

! Coefficient of the "Coriolis-type" term
if (icorio.eq.1) then
  ! Relative velocity formulation
  ccorio = 2.d0
else if (iturbo.eq.1) then
  ! Mixed relative/absolute velocity formulation
  ccorio = 1.d0
else
  ccorio = 0.d0
endif

!===============================================================================
! 2. Potential forces (pressure gradient and gravity)
!===============================================================================

!-------------------------------------------------------------------------------
! ---> Pressure gradient

! Allocate a work array for the gradient calculation
allocate(grad(3,ncelet))

iccocg = 1
inc    = 1

! For compressible flows, the new Pressure field is required
if (ippmod(icompf).ge.0) then
  iprev = 0
! For incompressible flows, keep the pressure at time n
! in case of PISO algorithm
else
  iprev = 1
endif

if (icavit.lt.0) then
  call field_gradient_potential(ivarfl(ipr), iprev, imrgra, inc,    &
                                iccocg, iphydr,                     &
                                frcxt, grad)

else

  ! Cavitating flows: consistency of the gradient with the diffusive flux scheme
  ! of the correction step

  call field_get_coefa_s (ivarfl(ipr), coefa_p)
  call field_get_coefb_s (ivarfl(ipr), coefb_p)

  allocate(xinvro(ncelet))

  do iel = 1, ncel
    xinvro(iel) = 1.d0/crom(iel)
  enddo

  iccocg = 1
  inc    = 1
  nswrgp = nswrgr(ipr)
  imligp = imligr(ipr)
  iwarnp = iwarni(ipr)
  epsrgp = epsrgr(ipr)
  climgp = climgr(ipr)
  extrap = extrag(ipr)

  call gradient_weighted_s(ivarfl(ipr), imrgra, inc, iccocg, nswrgp, imligp,  &
                           iwarnp, epsrgp, climgp, extrap,                    &
                           cvara_pr, xinvro, coefa_p, coefb_p,                &
                           grad )

  deallocate(xinvro)

endif

!    Calcul des efforts aux parois (partie 2/5), si demande
!    La pression a la face est calculee comme dans gradrc/gradmc
!    et on la transforme en pression totale
!    On se limite a la premiere iteration (pour faire simple par
!      rapport a la partie issue de condli, hors boucle)
if (iforbr.ge.0 .and. iterns.eq.1) then
  call field_get_coefa_s (ivarfl(ipr), coefa_p)
  call field_get_coefb_s (ivarfl(ipr), coefb_p)
  do ifac = 1, nfabor
    iel = ifabor(ifac)
    diipbx = diipb(1,ifac)
    diipby = diipb(2,ifac)
    diipbz = diipb(3,ifac)
    pip = cvara_pr(iel) &
        + diipbx*grad(1,iel) + diipby*grad(2,iel) + diipbz*grad(3,iel)
    pfac = coefa_p(ifac) +coefb_p(ifac)*pip
    pfac1= cvara_pr(iel)                                              &
         +(cdgfbo(1,ifac)-xyzcen(1,iel))*grad(1,iel)              &
         +(cdgfbo(2,ifac)-xyzcen(2,iel))*grad(2,iel)              &
         +(cdgfbo(3,ifac)-xyzcen(3,iel))*grad(3,iel)
    pfac = coefb_p(ifac)*(extrag(ipr)*pfac1                       &
         +(1.d0-extrag(ipr))*pfac)                                &
         +(1.d0-coefb_p(ifac))*pfac                               &
         + ro0*(gx*(cdgfbo(1,ifac)-xyzp0(1))                      &
         + gy*(cdgfbo(2,ifac)-xyzp0(2))                           &
         + gz*(cdgfbo(3,ifac)-xyzp0(3)) )                         &
         - pred0
    do isou = 1, 3
      forbr(isou,ifac) = forbr(isou,ifac) + pfac*surfbo(isou,ifac)
    enddo
  enddo
endif

!-------------------------------------------------------------------------------
! ---> RESIDU DE NORMALISATION POUR RESOLP
!     Test d'un residu de normalisation de l'etape de pression
!       plus comprehensible = div(rho u* + dt gradP^(n))-Gamma
!       i.e. second membre du systeme en pression hormis la partie
!       pression (sinon a convergence, on tend vers 0)
!       Represente les termes que la pression doit equilibrer
!     On calcule ici div(rho dt/rho gradP^(n)) et on complete a la fin
!       avec  div(rho u*)
!     Pour grad P^(n) on suppose que des CL de Neumann homogenes
!       s'appliquent partout : on peut donc utiliser les CL de la
!       vitesse pour u*+dt/rho gradP^(n). Comme on calcule en deux fois,
!       on utilise les CL de vitesse homogenes pour dt/rho gradP^(n)
!       ci-dessous et les CL de vitesse completes pour u* a la fin.

if (iappel.eq.1.and.irnpnw.eq.1) then

!     Calcul de dt/rho*grad P
  do iel = 1, ncel
    dtsrom = dt(iel)/crom(iel)
    trav(1,iel) = grad(1,iel)*dtsrom
    trav(2,iel) = grad(2,iel)*dtsrom
    trav(3,iel) = grad(3,iel)*dtsrom
  enddo

  if (irangp.ge.0.or.iperio.eq.1) then
    call synvin(trav)
  endif

!     Calcul de rho dt/rho*grad P.n aux faces
!       Pour gagner du temps, on ne reconstruit pas.
  itypfl = 1
  ! Cavitation algorithm: the pressure step corresponds to the
  ! correction of the volumetric flux, not the mass flux
  if (icavit.ge.0)  itypfl = 0
  init   = 1
  inc    = 0
  iflmb0 = 1
  nswrp  = 1
  imligp = imligr(iu )
  iwarnp = iwarni(ipr)
  epsrgp = epsrgr(iu )
  climgp = climgr(iu )

  call inimav                                                     &
 ( ivarfl(iu)      , itypfl ,                                     &
   iflmb0 , init   , inc    , imrgra , nswrp  , imligp ,          &
   iwarnp ,                                                       &
   epsrgp , climgp ,                                              &
   crom   , brom   ,                                              &
   trav   ,                                                       &
   coefav , coefbv ,                                              &
   viscf  , viscb  )

  ! Compute div(rho dt/rho*grad P)
  allocate(xnormp(ncelet))

  init = 1
  call divmas(init,viscf,viscb,xnormp)

!-- Volumic Gamma source term adding for volumic mass flow rate
  if (ncesmp.gt.0) then
    do ii = 1, ncesmp
      iel = icetsm(ii)
      xnormp(iel) = xnormp(iel) - cell_f_vol(iel)*smacel(ii,ipr)
    enddo
  endif

!-- Surface Gamma source term adding for surface condensation modelling
  if (nfbpcd.gt.0) then
    do ii = 1, nfbpcd
      ifac= ifbpcd(ii)
      iel = ifabor(ifac)
      xnormp(iel) = xnormp(iel) - surfbn(ifac) * spcond(ii,ipr)
    enddo
  endif

! --- volume Gamma source term adding for volume condensation modelling
  if (icond.eq.1) then
    allocate(surfbm(ncelet))
    surfbm(:) = 0.d0

    do ii = 1, ncmast
      iel= ltmast(ii)
      surfbm(iel) = s_metal*volume(iel)/voltot
      xnormp(iel) = xnormp(iel)  - surfbm(iel)*svcond(iel,ipr)
    enddo

    deallocate(surfbm)
  endif

  if (idilat.ge.4) then
    call field_get_val_s(iprpfl(iustdy(itsrho)), cpro_tsrho)
  endif

  ! Dilatable mass conservative algorithm
  if (idilat.eq.2) then
    do iel = 1, ncel
      drom = crom(iel) - croma(iel)
      xnormp(iel) = xnormp(iel) + drom*cell_f_vol(iel)/dt(iel)
    enddo
  ! Semi-analytic weakly compressible algorithm add + 1/rho Drho/Dt
  else if (idilat.eq.4)then
    do iel = 1, ncel
      xnormp(iel) = xnormp(iel) + cpro_tsrho(iel)/crom(iel)
    enddo
  else if (idilat.eq.5) then
    do iel = 1, ncel
      xnormp(iel) = xnormp(iel) + cpro_tsrho(iel)
    enddo
  endif

  ! Cavitation source term
  if (icavit.gt.0) then
    do iel = 1, ncel
      xnormp(iel) = xnormp(iel) -cell_f_vol(iel)*gamcav(iel)*(1.d0/rov - 1.d0/rol)
    enddo
  endif

!     On conserve XNORMP, on complete avec u* a la fin et
!       on le transfere a resopv

endif


!     Au premier appel, TRAV est construit directement ici.
!     Au second  appel (estimateurs), TRAV contient deja
!       l'increment temporel.
!     On pourrait fusionner en initialisant TRAV a zero
!       avant le premier appel, mais ca fait des operations en plus.

!     Remarques :
!       rho g sera a l'ordre 2 s'il est extrapole.
!       si on itere sur navsto, ca ne sert a rien de recalculer rho g a
!         chaque fois (ie on pourrait le passer dans trava) mais ce n'est
!         pas cher.
if (iappel.eq.1) then
  if (iphydr.eq.1) then
    do iel = 1, ncel
      trav(1,iel) = (frcxt(1 ,iel) - grad(1,iel)) * cell_f_vol(iel)
      trav(2,iel) = (frcxt(2 ,iel) - grad(2,iel)) * cell_f_vol(iel)
      trav(3,iel) = (frcxt(3 ,iel) - grad(3,iel)) * cell_f_vol(iel)
    enddo
  else if (iphydr.eq.2) then
    do iel = 1, ncel
      rom = crom(iel)
      trav(1,iel) = (rom*gx - grdphd(1,iel) - grad(1,iel)) * cell_f_vol(iel)
      trav(2,iel) = (rom*gy - grdphd(2,iel) - grad(2,iel)) * cell_f_vol(iel)
      trav(3,iel) = (rom*gz - grdphd(3,iel) - grad(3,iel)) * cell_f_vol(iel)
    enddo
  else if (ippmod(icompf).ge.0) then
    do iel = 1, ncel
      rom = crom(iel)
      trav(1,iel) = (rom*gx - grad(1,iel)) * cell_f_vol(iel)
      trav(2,iel) = (rom*gy - grad(2,iel)) * cell_f_vol(iel)
      trav(3,iel) = (rom*gz - grad(3,iel)) * cell_f_vol(iel)
    enddo
  else
    do iel = 1, ncel
      drom = (crom(iel)-ro0)
      trav(1,iel) = (drom*gx - grad(1,iel) ) * cell_f_vol(iel)
      trav(2,iel) = (drom*gy - grad(2,iel) ) * cell_f_vol(iel)
      trav(3,iel) = (drom*gz - grad(3,iel) ) * cell_f_vol(iel)
    enddo
  endif

else if(iappel.eq.2) then

  if (iphydr.eq.1) then
    do iel = 1, ncel
      trav(1,iel) = trav(1,iel) + (frcxt(1 ,iel) - grad(1,iel))*cell_f_vol(iel)
      trav(2,iel) = trav(2,iel) + (frcxt(2 ,iel) - grad(2,iel))*cell_f_vol(iel)
      trav(3,iel) = trav(3,iel) + (frcxt(3 ,iel) - grad(3,iel))*cell_f_vol(iel)
    enddo
  else if (iphydr.eq.2) then
    do iel = 1, ncel
      rom = crom(iel)
      trav(1,iel) = trav(1,iel) + (rom*gx - grdphd(1,iel) - grad(1,iel))*cell_f_vol(iel)
      trav(2,iel) = trav(2,iel) + (rom*gy - grdphd(2,iel) - grad(2,iel))*cell_f_vol(iel)
      trav(3,iel) = trav(3,iel) + (rom*gz - grdphd(3,iel) - grad(3,iel))*cell_f_vol(iel)
    enddo
  else
    do iel = 1, ncel
      drom = (crom(iel)-ro0)
      trav(1,iel) = trav(1,iel) + (drom*gx - grad(1,iel))*cell_f_vol(iel)
      trav(2,iel) = trav(2,iel) + (drom*gy - grad(2,iel))*cell_f_vol(iel)
      trav(3,iel) = trav(3,iel) + (drom*gz - grad(3,iel))*cell_f_vol(iel)
    enddo
  endif

endif

! Free memory
deallocate(grad)


!   Pour IAPPEL = 1 (ie appel standard sans les estimateurs)
!     TRAV rassemble les termes sources  qui seront recalcules
!       a toutes les iterations sur navsto
!     Si on n'itere pas sur navsto et qu'on n'extrapole pas les
!       termes sources, TRAV contient tous les termes sources
!       jusqu'au basculement dans SMBR
!     A ce niveau, TRAV contient -grad P et rho g
!       P est suppose pris a n+1/2
!       rho est eventuellement interpole a n+1/2


!-------------------------------------------------------------------------------
! ---> INITIALISATION DU TABLEAU TRAVA et terme source AU PREMIER PASSAGE
!     (A LA PREMIERE ITER SUR NAVSTO)

!     TRAVA rassemble les termes sources qu'il suffit de calculer
!       a la premiere iteration sur navsto quand il y a plusieurs iter.
!     Quand il n'y a qu'une iter, on cumule directement dans TRAV
!       ce qui serait autrement alle dans TRAVA
!     Les termes sources explicites serviront
!       pour le pas de temps suivant en cas d'extrapolation (plusieurs
!       iter sur navsto ou pas)

!     A la premiere iter sur navsto
if (iterns.eq.1) then

  ! Si on   extrapole     les T.S. : -theta*valeur precedente
  if (isno2t.gt.0) then
    ! S'il n'y a qu'une    iter : TRAV  incremente
    if (nterup.eq.1) then
      do iel = 1, ncel
        do ii = 1, ndim
          trav (ii,iel) = trav (ii,iel) - thets*c_st_vel(ii,iel)
        enddo
      enddo
      ! S'il   y a plusieurs iter : TRAVA initialise
    else
      do iel = 1, ncel
        do ii = 1, ndim
          trava(ii,iel) = - thets*c_st_vel(ii,iel)
        enddo
      enddo
    endif
    ! Et on initialise le terme source pour le remplir ensuite
    do iel = 1, ncel
      do ii = 1, ndim
        c_st_vel(ii,iel) = 0.d0
      enddo
    enddo

  ! Si on n'extrapole pas les T.S.
  else
    ! S'il   y a plusieurs iter : TRAVA initialise
    !  sinon TRAVA n'existe pas
    if(nterup.gt.1) then
      do ii = 1, ndim
        do iel = 1, ncel
          trava(ii,iel)  = 0.d0
        enddo
      enddo
    endif
  endif

endif

!-------------------------------------------------------------------------------
! Initialization of the implicit terms

if (iappel.eq.1) then

  ! Low Mach compressible Algos
  if (idilat.gt.1.or.ippmod(icompf).ge.0) then
    call field_get_val_prev_s(icrom, pcrom)

  ! Cavitation
  else if (icavit.ge.0) then
    call field_get_val_s(icroaa, pcrom)

  ! Standard algo
  else

    call field_get_val_s(icrom, pcrom)
  endif

  do iel = 1, ncel
    do isou = 1, 3
      fimp(isou,isou,iel) = istat(iu)*pcrom(iel)/dt(iel)*cell_f_vol(iel)
      do jsou = 1, 3
        if(jsou.ne.isou) fimp(isou,jsou,iel) = 0.d0
      enddo
    enddo
  enddo

!     Le remplissage de FIMP est toujours indispensable,
!       meme si on peut se contenter de n'importe quoi pour IAPPEL=2.
else
  do iel = 1, ncel
    do isou = 1, 3
      do jsou = 1, 3
        fimp(isou,jsou,iel) = 0.d0
      enddo
    enddo
  enddo
endif

!-------------------------------------------------------------------------------
! ---> 2/3 RHO * GRADIENT DE K SI k-epsilon ou k-omega
!      NB : ON NE PREND PAS LE GRADIENT DE (RHO K), MAIS
!           CA COMPLIQUERAIT LA GESTION DES CL ...
!     On peut se demander si l'extrapolation en temps sert a
!       quelquechose

!     Ce terme explicite est calcule une seule fois,
!       a la premiere iter sur navsto : il est stocke dans un champ si on
!       doit l'extrapoler en temps ; il va dans TRAVA si on n'extrapole
!       pas mais qu'on itere sur navsto. Il va dans TRAV si on
!       n'extrapole pas et qu'on n'itere pas sur navsto.
if(     (itytur.eq.2 .or. itytur.eq.5 .or. iturb.eq.60) &
   .and. igrhok.eq.1 .and. iterns.eq.1) then

  ! Allocate a work array for the gradient calculation
  allocate(grad(3,ncelet))

  iccocg = 1
  iprev  = 1
  inc    = 1

  call field_gradient_scalar(ivarfl(ik), iprev, imrgra, inc,      &
                             iccocg,                              &
                             grad)

  d2s3 = 2.d0/3.d0

  ! Si on extrapole les termes source en temps
  if (isno2t.gt.0) then
    ! Calcul de rho^n grad k^n      si rho non extrapole
    !           rho^n grad k^n      si rho     extrapole

    call field_get_val_s(icrom, crom)
    call field_get_val_prev_s(icrom, croma)
    do iel = 1, ncel
      romvom = -croma(iel)*cell_f_vol(iel)*d2s3
      do isou = 1, 3
        c_st_vel(isou,iel) = c_st_vel(isou,iel)+grad(isou,iel)*romvom
      enddo
    enddo
  ! Si on n'extrapole pas les termes sources en temps : TRAV ou TRAVA
  else
    if(nterup.eq.1) then
      do iel = 1, ncel
        romvom = -crom(iel)*cell_f_vol(iel)*d2s3
        do isou = 1, 3
          trav(isou,iel) = trav(isou,iel) + grad(isou,iel) * romvom
        enddo
      enddo
    else
      do iel = 1, ncel
        romvom = -crom(iel)*cell_f_vol(iel)*d2s3
        do isou = 1, 3
          trava(isou,iel) = trava(isou,iel) + grad(isou,iel) * romvom
        enddo
      enddo
    endif
  endif

  ! Calcul des efforts aux parois (partie 3/5), si demande
  if (iforbr.ge.0) then
    call field_get_coefa_s (ivarfl(ik), coefa_k)
    call field_get_coefb_s (ivarfl(ik), coefb_k)
    do ifac = 1, nfabor
      iel = ifabor(ifac)
      diipbx = diipb(1,ifac)
      diipby = diipb(2,ifac)
      diipbz = diipb(3,ifac)
      xkb = cvara_k(iel) + diipbx*grad(1,iel)                      &
           + diipby*grad(2,iel) + diipbz*grad(3,iel)
      xkb = coefa_k(ifac)+coefb_k(ifac)*xkb
      xkb = d2s3*crom(iel)*xkb
      do isou = 1, 3
        forbr(isou,ifac) = forbr(isou,ifac) + xkb*surfbo(isou,ifac)
      enddo
    enddo
  endif

  ! Free memory
  deallocate(grad)

endif


!-------------------------------------------------------------------------------
! ---> Transpose of velocity gradient in the diffusion term

!     These terms are taken into account in bilscv.
!     We only compute here the secondary viscosity.

if (ivisse.eq.1) then

  call visecv(secvif, secvib)

endif

!-------------------------------------------------------------------------------
! ---> Head losses
!      (if iphydr=1 this term has already been taken into account)

! ---> Explicit part
if ((ncepdp.gt.0).and.(iphydr.ne.1)) then

  ! Les termes diagonaux sont places dans TRAV ou TRAVA,
  !   La prise en compte de uvwk a partir de la seconde iteration
  !   est faite directement dans coditv.
  if (iterns.eq.1) then

    ! On utilise temporairement TRAV comme tableau de travail.
    ! Son contenu est stocke dans W7, W8 et W9 jusqu'apres tspdcv
    do iel = 1,ncel
      w7(iel) = trav(1,iel)
      w8(iel) = trav(2,iel)
      w9(iel) = trav(3,iel)
      trav(1,iel) = 0.d0
      trav(2,iel) = 0.d0
      trav(3,iel) = 0.d0
    enddo

    call tspdcv(ncepdp, icepdc, vela, ckupdc, trav)

    ! Si on itere sur navsto, on utilise TRAVA ; sinon TRAV
    if(nterup.gt.1) then
      do iel = 1, ncel
        trava(1,iel) = trava(1,iel) + trav(1,iel)
        trava(2,iel) = trava(2,iel) + trav(2,iel)
        trava(3,iel) = trava(3,iel) + trav(3,iel)
        trav(1,iel)  = w7(iel)
        trav(2,iel)  = w8(iel)
        trav(3,iel)  = w9(iel)
      enddo
    else
      do iel = 1, ncel
        trav(1,iel)  = w7(iel) + trav(1,iel)
        trav(2,iel)  = w8(iel) + trav(2,iel)
        trav(3,iel)  = w9(iel) + trav(3,iel)
      enddo
    endif
  endif

endif

! ---> Implicit part

!  At the second call, fimp is not needed anymore
if (iappel.eq.1) then
  if (ncepdp.gt.0) then
    ! The theta-scheme for the head loss is the same as the other terms
    thetap = thetav(iu)
    do ielpdc = 1, ncepdp
      iel = icepdc(ielpdc)
      romvom = crom(iel)*cell_f_vol(iel)*thetap

      ! Diagonal part
      do isou = 1, 3
        fimp(isou,isou,iel) = fimp(isou,isou,iel) + romvom*ckupdc(ielpdc,isou)
      enddo
      ! Extra-diagonal part
      cpdc12 = ckupdc(ielpdc,4)
      cpdc23 = ckupdc(ielpdc,5)
      cpdc13 = ckupdc(ielpdc,6)

      fimp(1,2,iel) = fimp(1,2,iel) + romvom*cpdc12
      fimp(2,1,iel) = fimp(2,1,iel) + romvom*cpdc12
      fimp(1,3,iel) = fimp(1,3,iel) + romvom*cpdc13
      fimp(3,1,iel) = fimp(3,1,iel) + romvom*cpdc13
      fimp(2,3,iel) = fimp(2,3,iel) + romvom*cpdc23
      fimp(3,2,iel) = fimp(3,2,iel) + romvom*cpdc23
    enddo
  endif
endif


!-------------------------------------------------------------------------------
! ---> Coriolis force
!     (if iphydr=1 then this term is already taken into account)

! --->  Explicit part

if ((icorio.eq.1.or.iturbo.eq.1) .and. iphydr.ne.1) then

  ! A la premiere iter sur navsto, on ajoute la partie issue des
  ! termes explicites
  if (iterns.eq.1) then

    ! Si on n'itere pas sur navsto : TRAV
    if (nterup.eq.1) then

      call field_get_val_s(icrom, crom)

      do iel = 1, ncel
        romvom = -ccorio*crom(iel)*cell_f_vol(iel)
        call add_coriolis_v(irotce(iel), romvom, vela(:,iel), trav(:,iel))
      enddo

    ! Si on itere sur navsto : TRAVA
    else

      do iel = 1, ncel
        romvom = -ccorio*crom(iel)*cell_f_vol(iel)
        call add_coriolis_v(irotce(iel), romvom, vela(:,iel), trava(:,iel))
      enddo

    endif
  endif
endif

! --->  Implicit part

!  At the second call, fimp is not needed anymore
if(iappel.eq.1) then
  if (icorio.eq.1 .or. iturbo.eq.1) then
    ! The theta-scheme for the Coriolis term is the same as the other terms
    thetap = thetav(iu)

    do iel = 1, ncel
      romvom = -ccorio*crom(iel)*cell_f_vol(iel)*thetap
      call add_coriolis_t(irotce(iel), romvom, fimp(:,:,iel))
    enddo

  endif
endif

!-------------------------------------------------------------------------------
! ---> - Divergence of tensor Rij

if(itytur.eq.3.and.iterns.eq.1) then

  allocate(rij(6,ncelet))
  if(irijco.eq.1) then !TODO change index of rij
    do iel = 1, ncelet
      rij(1,iel) = cvara_rij(1,iel)
      rij(2,iel) = cvara_rij(2,iel)
      rij(3,iel) = cvara_rij(3,iel)
      rij(4,iel) = cvara_rij(4,iel)
      rij(5,iel) = cvara_rij(5,iel)
      rij(6,iel) = cvara_rij(6,iel)
    enddo
  else
    do iel = 1, ncelet
      rij(1,iel) = cvara_r11(iel)
      rij(2,iel) = cvara_r22(iel)
      rij(3,iel) = cvara_r33(iel)
      rij(4,iel) = cvara_r12(iel)
      rij(5,iel) = cvara_r23(iel)
      rij(6,iel) = cvara_r13(iel)
    enddo
  endif
! --- Boundary conditions on the components of the tensor Rij

  allocate(coefat(6,nfabor))
  if(irijco.eq.1) then
    call field_get_coefad_v(ivarfl(irij),coefap)
    coefat = coefap
  else
    call field_get_coefad_s(ivarfl(ir11),coef1)
    call field_get_coefad_s(ivarfl(ir22),coef2)
    call field_get_coefad_s(ivarfl(ir33),coef3)
    call field_get_coefad_s(ivarfl(ir12),coef4)
    call field_get_coefad_s(ivarfl(ir23),coef5)
    call field_get_coefad_s(ivarfl(ir13),coef6)
    do ifac = 1, nfabor
      coefat(1,ifac) = coef1(ifac)
      coefat(2,ifac) = coef2(ifac)
      coefat(3,ifac) = coef3(ifac)
      coefat(4,ifac) = coef4(ifac)
      coefat(5,ifac) = coef5(ifac)
      coefat(6,ifac) = coef6(ifac)
    enddo
  endif

  allocate(coefbt(6,6,nfabor))
  do ifac = 1, nfabor
    do ii = 1, 6
      do jj = 1, 6
        coefbt(jj,ii,ifac) = 0.d0
      enddo
    enddo
  enddo

  if(irijco.eq.1) then
    call field_get_coefbd_v(ivarfl(irij),coefbp)
    coefbt = coefbp
  else
    call field_get_coefbd_s(ivarfl(ir11),coef1)
    call field_get_coefbd_s(ivarfl(ir22),coef2)
    call field_get_coefbd_s(ivarfl(ir33),coef3)
    call field_get_coefbd_s(ivarfl(ir12),coef4)
    call field_get_coefbd_s(ivarfl(ir23),coef5)
    call field_get_coefbd_s(ivarfl(ir13),coef6)
    do ifac = 1, nfabor
      coefbt(1,1,ifac) = coef1(ifac)
      coefbt(2,2,ifac) = coef2(ifac)
      coefbt(3,3,ifac) = coef3(ifac)
      coefbt(4,4,ifac) = coef4(ifac)
      coefbt(5,5,ifac) = coef5(ifac)
      coefbt(6,6,ifac) = coef6(ifac)
    enddo
  endif

  ! Flux computation options
  f_id = -1
  init = 1;
  inc  = 1;
  iflmb0 = 0;
  nswrgp = nswrgr(ir11);
  imligp = imligr(ir11);
  iwarnp = iwarni(ir11);
  epsrgp = epsrgr(ir11);
  climgp = climgr(ir11);
  itypfl = 1;

  allocate(tflmas(3,nfac))
  allocate(tflmab(3,nfabor))

  call divrij &
 ( f_id   , itypfl ,                                              &
   iflmb0 , init   , inc    , imrgra , nswrgp , imligp ,          &
   iwarnp ,                                                       &
   epsrgp , climgp ,                                              &
   crom   , brom   ,                                              &
   rij    ,                                                       &
   coefat , coefbt ,                                              &
   tflmas , tflmab )

  deallocate(rij)
  deallocate(coefat, coefbt)

  !     Calcul des efforts aux bords (partie 5/5), si necessaire

  if (iforbr.ge.0) then
    do ifac = 1, nfabor
      do isou = 1, 3
        forbr(isou,ifac) = forbr(isou,ifac) + tflmab(isou,ifac)
      enddo
    enddo
  endif

  allocate(divt(3,ncelet))
  init = 1
  call divmat(init,tflmas,tflmab,divt)

  deallocate(tflmas, tflmab)

  ! (if iphydr=1 then this term is already taken into account)
  if (iphydr.ne.1.or.igprij.ne.1) then

    ! If extrapolation of source terms
    if (isno2t.gt.0) then
      do iel = 1, ncel
        do isou = 1, 3
          c_st_vel(isou,iel) = c_st_vel(isou,iel) - divt(isou,iel)
        enddo
      enddo

    ! No extrapolation of source terms
    else

      ! No PISO iteration
      if (nterup.eq.1) then
        do iel = 1, ncel
          do isou = 1, 3
            trav(isou,iel) = trav(isou,iel) - divt(isou,iel)
          enddo
        enddo
      ! PISO iterations
      else
        do iel = 1, ncel
          do isou = 1, 3
            trava(isou,iel) = trava(isou,iel) - divt(isou,iel)
          enddo
        enddo
      endif
    endif
  endif

endif


!-------------------------------------------------------------------------------
! ---> Face diffusivity for the velocity

if (idiff(iu).ge. 1) then

  call field_get_val_s(iprpfl(iviscl), viscl)
  call field_get_val_s(iprpfl(ivisct), visct)

  if (itytur.eq.3) then
    do iel = 1, ncel
      w1(iel) = viscl(iel)
    enddo
  else
    do iel = 1, ncel
      w1(iel) = viscl(iel) + idifft(iu)*visct(iel)
    enddo
  endif

  ! Scalar diffusivity (Default)
  if (idften(iu).eq.1) then

    call viscfa &
   ( imvisf ,                                                       &
     w1     ,                                                       &
     viscf  , viscb  )

    ! When using Rij-epsilon model with the option irijnu=1, the face
    ! viscosity for the Matrix (viscfi and viscbi) is increased
    if(itytur.eq.3.and.irijnu.eq.1) then

      do iel = 1, ncel
        w1(iel) = viscl(iel) + idifft(iu)*visct(iel)
      enddo

      call viscfa &
   ( imvisf ,                                                       &
     w1     ,                                                       &
     viscfi , viscbi )
    endif

  ! Tensorial diffusion of the velocity (in case of tensorial porosity)
  else if (idften(iu).eq.6) then

    do iel = 1, ncel
      do isou = 1, 3
        viscce(isou, iel) = w1(iel)
      enddo
      do isou = 4, 6
        viscce(isou, iel) = 0.d0
      enddo
    enddo

    call vistnv &
     ( imvisf ,                                                       &
       viscce ,                                                       &
       viscf  , viscb  )

    ! When using Rij-epsilon model with the option irijnu=1, the face
    ! viscosity for the Matrix (viscfi and viscbi) is increased
    if(itytur.eq.3.and.irijnu.eq.1) then

      do iel = 1, ncel
        w1(iel) = viscl(iel) + idifft(iu)*visct(iel)
      enddo

      do iel = 1, ncel
        do isou = 1, 3
          viscce(isou, iel) = w1(iel)
        enddo
        do isou = 4, 6
          viscce(isou, iel) = 0.d0
        enddo
      enddo

      call vistnv &
       ( imvisf ,                                                       &
         viscce ,                                                       &
         viscfi , viscbi )

    endif
  endif

! --- If no diffusion, viscosity is set to 0.
else

  do ifac = 1, nfac
    viscf(ifac) = 0.d0
  enddo
  do ifac = 1, nfabor
    viscb(ifac) = 0.d0
  enddo

  if(itytur.eq.3.and.irijnu.eq.1) then
    do ifac = 1, nfac
      viscfi(ifac) = 0.d0
    enddo
    do ifac = 1, nfabor
      viscbi(ifac) = 0.d0
    enddo
  endif

endif

!-------------------------------------------------------------------------------
! ---> Take external forces partially equilibrated with the pressure gradient
!      into account (only for the first call, the second one is dedicated
!      to error estimators)

if (iappel.eq.1.and.iphydr.eq.1.and.iterns.eq.1) then

! force ext au pas de temps precedent :
!     FRCXT a ete initialise a zero
!     (est deja utilise dans typecl, et est mis a jour a la fin
!     de navsto)

  do iel = 1, ncel

    ! External force variation between time step n and n+1
    ! (used in the correction step)
    drom = (crom(iel)-ro0)
    dfrcxt(1, iel) = drom*gx - frcxt(1, iel)
    dfrcxt(2, iel) = drom*gy - frcxt(2, iel)
    dfrcxt(3, iel) = drom*gz - frcxt(3, iel)
  enddo

  ! Add head losses
  if (ncepdp.gt.0) then
    do ielpdc = 1, ncepdp
      iel=icepdc(ielpdc)
      vit1   = vela(1,iel)
      vit2   = vela(2,iel)
      vit3   = vela(3,iel)
      cpdc11 = ckupdc(ielpdc,1)
      cpdc22 = ckupdc(ielpdc,2)
      cpdc33 = ckupdc(ielpdc,3)
      cpdc12 = ckupdc(ielpdc,4)
      cpdc23 = ckupdc(ielpdc,5)
      cpdc13 = ckupdc(ielpdc,6)
      dfrcxt(1 ,iel) = dfrcxt(1 ,iel) &
                    - crom(iel)*(cpdc11*vit1+cpdc12*vit2+cpdc13*vit3)
      dfrcxt(2 ,iel) = dfrcxt(2 ,iel) &
                    - crom(iel)*(cpdc12*vit1+cpdc22*vit2+cpdc23*vit3)
      dfrcxt(3 ,iel) = dfrcxt(3 ,iel) &
                    - crom(iel)*(cpdc13*vit1+cpdc23*vit2+cpdc33*vit3)
    enddo
  endif

  ! Add Coriolis force
  if (icorio.eq.1 .or. iturbo.eq.1) then
    do iel = 1, ncel
      rom = -ccorio*crom(iel)
      call add_coriolis_v(irotce(iel), rom, vela(:,iel), dfrcxt(:,iel))
    enddo
  endif

  ! Add -div( rho R) as external force
  if (itytur.eq.3.and.igprij.eq.1) then
    do iel = 1, ncel
      dvol = 1.d0/cell_f_vol(iel)
      do isou = 1, 3
        dfrcxt(isou, iel) = dfrcxt(isou, iel) - divt(isou, iel)*dvol
      enddo
    enddo
  endif


  if (irangp.ge.0.or.iperio.eq.1) then
    call synvin(dfrcxt)
  endif

endif

!===============================================================================
! 3. Solving of the 3x3xNcel coupled system
!===============================================================================


! ---> AU PREMIER APPEL,
!      MISE A ZERO DE L'ESTIMATEUR POUR LA VITESSE PREDITE
!      S'IL DOIT ETRE CALCULE

if (iappel.eq.1) then
  if (iestim(iespre).ge.0) then
    call field_get_val_s(iestim(iespre), c_estim)
    do iel = 1, ncel
      c_estim(iel) =  0.d0
    enddo
  endif
endif

! ---> AU DEUXIEME APPEL,
!      MISE A ZERO DE L'ESTIMATEUR TOTAL POUR NAVIER-STOKES
!      (SI ON FAIT UN DEUXIEME APPEL, ALORS IL DOIT ETRE CALCULE)

if (iappel.eq.2) then
  call field_get_val_s(iestim(iestot), c_estim)
  do iel = 1, ncel
    c_estim(iel) =  0.d0
  enddo
endif

!-------------------------------------------------------------------------------
! ---> User source terms

do iel = 1, ncel
  do isou = 1, 3
    tsexp(isou,iel) = 0.d0
    do jsou = 1, 3
      tsimp(isou,jsou,iel) = 0.d0
    enddo
  enddo
enddo

! The computation of esplicit and implicit source terms is performed
! at the first iteration only.
if (iterns.eq.1) then

  if (iihmpr.eq.1) then
    call uitsnv (vel, tsexp, tsimp)
  endif

  n_fans = cs_fan_n_fans()
  if (n_fans .gt. 0) then
    if (ntcabs.eq.ntpabs+1) then
      call debvtl(flumas, flumab, crom, brom)
    endif
    call tsvvtl(tsexp)
  endif

  call ustsnv &
 ( nvar   , nscal  , ncepdp , ncesmp ,                            &
   iu   ,                                                         &
   icepdc , icetsm , itypsm ,                                     &
   dt     ,                                                       &
   ckupdc , smacel , tsexp  , tsimp  )

  if (ibdtso(iu).gt.1.and.ntcabs.gt.ntinit &
      .and.(idtvar.eq.0.or.idtvar.eq.1)) then
    ! TODO: remove test on ntcabs and implemente a "proper" condition for
    ! initialization.
    f_id = ivarfl(iu)
    call cs_backward_differentiation_in_time(f_id, tsexp, tsimp)
  endif
  ! Skip first time step after restart if previous values have not been read.
  if (ibdtso(iu).lt.0) ibdtso(iu) = iabs(ibdtso(iu))

  ! Coupling between two Code_Saturne
  if (nbrcpl.gt.0) then
    !vectorial interleaved exchange
    call csccel(iu, vela, coefav, coefbv, tsexp)
  endif

  if (iphydr.eq.1.and.igpust.eq.1) then

    do iel = 1, ncel
      !FIXME when using porosity
      dvol = 1.d0/cell_f_vol(iel)
      do isou = 1, 3
        dfrcxt(isou, iel) = dfrcxt(isou, iel) + tsexp(isou, iel)*dvol
      enddo
    enddo

    if (irangp.ge.0.or.iperio.eq.1) then
      call synvin(dfrcxt)
    endif
  endif

endif

! if PISO sweeps are expected, implicit user sources terms are stored in ximpa
if (iterns.eq.1.and.nterup.gt.1) then
  do iel = 1, ncel
    do isou = 1, 3
      do jsou = 1, 3
        ximpa(isou,jsou,iel) = tsimp(isou,jsou,iel)
      enddo
    enddo
  enddo
endif

! ---> Explicit contribution due to implicit terms

if (iterns.eq.1) then
  if (nterup.gt.1) then
    do iel = 1, ncel
      do isou = 1, 3
        do jsou = 1, 3
          trava(isou,iel) = trava(isou,iel)                                  &
                          + tsimp(isou,jsou,iel)*vela(jsou,iel)
        enddo
      enddo
    enddo
  else
    do iel = 1, ncel
      do isou = 1, 3
        do jsou = 1, 3
          trav(isou,iel) = trav(isou,iel)                                    &
                         + tsimp(isou,jsou,iel)*vela(jsou,iel)
        enddo
      enddo
    enddo
  endif
endif

! At the first PISO iteration, explicit source terms are added
if (iterns.eq.1.and.(iphydr.ne.1.or.igpust.ne.1)) then
  ! If source terms are time-extrapolated, they are stored in fields
  if (isno2t.gt.0) then
    do iel = 1, ncel
      do isou = 1, 3
        c_st_vel(isou,iel) = c_st_vel(isou,iel) + tsexp(isou,iel)
      enddo
    enddo

  else
    ! If no PISO sweep
    if (nterup.eq.1) then
      do iel = 1, ncel
        do isou = 1, 3
          trav(isou,iel) = trav(isou,iel) + tsexp(isou,iel)
        enddo
      enddo
    ! If PISO sweeps
    else
      do iel = 1, ncel
        do isou = 1, 3
          trava(isou,iel) = trava(isou,iel) + tsexp(isou,iel)
        enddo
      enddo
    endif
  endif
endif

! ---> Implicit terms
if (iappel.eq.1) then
  ! If source terms are time-extrapolated
  if (isno2t.gt.0) then
    thetap = thetav(iu)
    if (iterns.gt.1) then
      do iel = 1, ncel
        do isou = 1, 3
          do jsou = 1, 3
            fimp(isou,jsou,iel) = fimp(isou,jsou,iel)                      &
                                - ximpa(isou,jsou,iel)*thetap
          enddo
        enddo
      enddo
    else
      do iel = 1, ncel
        do isou = 1, 3
          do jsou = 1, 3
            fimp(isou,jsou,iel) = fimp(isou,jsou,iel)                      &
                                - tsimp(isou,jsou,iel)*thetap
          enddo
        enddo
      enddo
    endif
  else
    if (iterns.gt.1) then
      do iel = 1, ncel
        do isou = 1, 3
          do jsou = 1, 3
            fimp(isou,jsou,iel) = fimp(isou,jsou,iel)                      &
                                + max(-ximpa(isou,jsou,iel),zero)
          enddo
        enddo
      enddo
    else
      do iel = 1, ncel
        do isou = 1, 3
          do jsou = 1, 3
            fimp(isou,jsou,iel) = fimp(isou,jsou,iel)                      &
                                + max(-tsimp(isou,jsou,iel),zero)
          enddo
        enddo
      enddo
    endif
  endif
endif

!-------------------------------------------------------------------------------
! --->  Mass source terms

if (ncesmp.gt.0) then

!     On calcule les termes Gamma (uinj - u)
!       -Gamma u a la premiere iteration est mis dans
!          TRAV ou TRAVA selon qu'on itere ou non sur navsto
!       Gamma uinj a la premiere iteration est placee dans W1
!       ROVSDT a chaque iteration recoit Gamma
  allocate(gavinj(3,ncelet))
  if (nterup.eq.1) then
    call catsmv &
  ( ncelet , ncel , ncesmp , iterns , isno2t,                   &
    icetsm , itypsm(1,iu),                                      &
    cell_f_vol    , vela , smacel(1,iu) , smacel(1,ipr) ,       &
    trav   , fimp , gavinj )
  else
    call catsmv &
  ( ncelet , ncel , ncesmp , iterns , isno2t,                   &
    icetsm , itypsm(1,iu),                                      &
    cell_f_vol    , vela , smacel(1,iu) , smacel(1,ipr) ,       &
    trava  , fimp  , gavinj )
  endif

  ! At the first PISO iteration, the explicit part "Gamma u^{in}" is added
  if (iterns.eq.1) then
    ! If source terms are extrapolated, stored in fields
    if(isno2t.gt.0) then
      do iel = 1, ncel
        do isou = 1, 3
          c_st_vel(isou,iel) = c_st_vel(isou,iel) + gavinj(isou,iel)
        enddo
      enddo

    else
      ! If no PISO iteration: in trav
      if (nterup.eq.1) then
        do iel = 1,ncel
          do isou = 1, 3
            trav(isou,iel)  = trav(isou,iel) + gavinj(isou,iel)
          enddo
        enddo
      ! Otherwise, in trava
      else
        do iel = 1,ncel
          do isou = 1, 3
            trava(isou,iel) = trava(isou,iel) + gavinj(isou,iel)
          enddo
        enddo
      endif
    endif
  endif

  deallocate(gavinj)

endif

! ---> Right Han Side initialization

! If source terms are extrapolated in time
if (isno2t.gt.0) then
  thetp1 = 1.d0 + thets
  ! If no PISO iteration: trav
  if (nterup.eq.1) then
    do iel = 1, ncel
      do isou = 1, 3
        smbr(isou,iel) = trav(isou,iel) + thetp1*c_st_vel(isou,iel)
      enddo
    enddo

  else
    do iel = 1, ncel
      do isou = 1, 3
        smbr(isou,iel) = trav(isou,iel) + trava(isou,iel)       &
                       + thetp1*c_st_vel(isou,iel)
      enddo
    enddo
  endif

! No time extrapolation
else
  ! No PISO iteration
  if (nterup.eq.1) then
    do iel = 1, ncel
      do isou = 1, 3
        smbr(isou,iel) = trav(isou,iel)
      enddo
    enddo
  ! PISO iterations
  else
    do iel = 1, ncel
      do isou = 1, 3
        smbr(isou,iel) = trav(isou,iel) + trava(isou,iel)
      enddo
    enddo
  endif
endif

! ---> LAGRANGIEN : COUPLAGE RETOUR

!     L'ordre 2 sur les termes issus du lagrangien necessiterait de
!       decomposer TSLAGR(IEL,ISOU) en partie implicite et
!       explicite, comme c'est fait dans ustsnv.
!     Pour le moment, on n'y touche pas.

if (iilagr.eq.2 .and. ltsdyn.eq.1)  then

  do iel = 1, ncel
    do isou = 1, 3
      smbr(isou,iel) = smbr(isou,iel) + tslagr(iel,itsvx+isou-1)
    enddo
  enddo
  ! At the second call, fimp is unused
  if(iappel.eq.1) then
    do iel = 1, ncel
      do isou = 1, 3
        fimp(isou,isou,iel) = fimp(isou,isou,iel) + max(-tslagr(iel,itsli),zero)
      enddo
    enddo
  endif

endif

! ---> Electric Arc (Laplace Force)
!      (No 2nd order in time yet)
if (ippmod(ielarc).ge.1) then
  call field_get_val_v_by_name('laplace_force', lapla)
  do iel = 1, ncel
      smbr(1,iel) = smbr(1,iel) + cell_f_vol(iel) * lapla(1,iel)
      smbr(2,iel) = smbr(2,iel) + cell_f_vol(iel) * lapla(2,iel)
      smbr(3,iel) = smbr(3,iel) + cell_f_vol(iel) * lapla(3,iel)
  enddo
endif

! Solver parameters
iconvp = iconv (iu)
idiffp = idiff (iu)
ndircp = ndircl(iu)
nswrsp = nswrsm(iu)
nswrgp = nswrgr(iu)
imligp = imligr(iu)
ircflp = ircflu(iu)
ischcp = ischcv(iu)
isstpp = isstpc(iu)
idftnp = idften(iu)
iswdyp = iswdyn(iu)
iwarnp = iwarni(iu)
blencp = blencv(iu)
epsilp = epsilo(iu)
epsrsp = epsrsm(iu)
epsrgp = epsrgr(iu)
climgp = climgr(iu)
extrap = extrag(iu)
relaxp = relaxv(iu)
thetap = thetav(iu)

if (ippmod(icompf).ge.0) then
  ! impose boundary convective flux at some faces (face indicator icvfli)
  icvflb = 1
else
  ! all boundary convective flux with upwind
  icvflb = 0
endif

if (iappel.eq.1) then

  iescap = iescal(iespre)

  if (iterns.eq.1) then

    ! Warning: in case of convergence estimators, eswork give the estimator
    ! of the predicted velocity
    call coditv &
 ( idtvar , iu     , iconvp , idiffp , ndircp ,                   &
   imrgra , nswrsp , nswrgp , imligp , ircflp , ivisse ,          &
   ischcp , isstpp , iescap , idftnp , iswdyp ,                   &
   iwarnp ,                                                       &
   blencp , epsilp , epsrsp , epsrgp , climgp ,                   &
   relaxp , thetap ,                                              &
   vela   , vela   ,                                              &
   coefav , coefbv , cofafv , cofbfv ,                            &
   flumas , flumab ,                                              &
   viscfi , viscbi , viscf  , viscb  , secvif , secvib ,          &
   icvflb , icvfli ,                                              &
   fimp   ,                                                       &
   smbr   ,                                                       &
   vel    ,                                                       &
   eswork )

  else if(iterns.gt.1) then

    call coditv &
 ( idtvar , iu     , iconvp , idiffp , ndircp ,                   &
   imrgra , nswrsp , nswrgp , imligp , ircflp , ivisse ,          &
   ischcp , isstpp , iescap , idftnp , iswdyp ,                   &
   iwarnp ,                                                       &
   blencp , epsilp , epsrsp , epsrgp , climgp ,                   &
   relaxp , thetap ,                                              &
   vela   , uvwk   ,                                              &
   coefav , coefbv , cofafv , cofbfv ,                            &
   flumas , flumab ,                                              &
   viscfi , viscbi , viscf  , viscb  , secvif , secvib ,          &
   icvflb , icvfli ,                                              &
   fimp   ,                                                       &
   smbr   ,                                                       &
   vel    ,                                                       &
   eswork )

  endif

  ! Velocity-pression coupling: compute the vector T, stored in tpucou,
  !  coditv is called, only one sweep is done, and tpucou is initialized
  !  by 0. so that the advection/diffusion added by bilscv is 0.
  !  nswrsp = -1 indicated that only one sweep is required and inc=0
  !  for boundary contitions on the weight matrix.
  if (ipucou.eq.1) then

    ! Allocate temporary arrays for the velocity-pressure resolution
    allocate(vect(3,ncelet))

    nswrsp = -1
    do iel = 1, ncel
      do isou = 1, 3
        smbr(isou,iel) = cell_f_vol(iel)
      enddo
    enddo
    do iel = 1, ncelet
      do isou = 1, 3
        vect(isou,iel) = 0.d0
      enddo
    enddo
    iescap = 0

    ! We do not take into account transpose of grad
    ivisep = 0

    call coditv &
 ( idtvar , iu     , iconvp , idiffp , ndircp ,                   &
   imrgra , nswrsp , nswrgp , imligp , ircflp , ivisep ,          &
   ischcp , isstpp , iescap , idftnp , iswdyp ,                   &
   iwarnp ,                                                       &
   blencp , epsilp , epsrsp , epsrgp , climgp ,                   &
   relaxp , thetap ,                                              &
   vect   , vect   ,                                              &
   coefav , coefbv , cofafv , cofbfv ,                            &
   flumas , flumab ,                                              &
   viscfi , viscbi , viscf  , viscb  , secvif , secvib ,          &
   icvflb , ivoid  ,                                              &
   fimp   ,                                                       &
   smbr   ,                                                       &
   vect   ,                                                       &
   rvoid  )

    do iel = 1, ncelet
      rom = crom(iel)
      do isou = 1, 3
        tpucou(isou,iel) = rom*vect(isou,iel)
      enddo
      do isou = 4, 6
        tpucou(isou,iel) = 0.d0
      enddo
    enddo

    ! Free memory
    deallocate(vect)

  endif

  ! ---> The estimator on the predicted velocity is summed up over the components
  if (iestim(iespre).ge.0) then
    call field_get_val_s(iestim(iespre), c_estim)
    do iel = 1, ncel
      do isou = 1, 3
        c_estim(iel) =  c_estim(iel) + eswork(isou,iel)
      enddo
    enddo
  endif


! ---> End of the construction of the total estimator:
!       RHS resiudal of (U^{n+1}, P^{n+1}) + rho*volume*(U^{n+1} - U^n)/dt
else if (iappel.eq.2) then

  inc = 1
  ! Pas de relaxation en stationnaire
  idtva0 = 0
  imasac = 0

  call bilscv &
 ( idtva0 , iu     , iconvp , idiffp , nswrgp , imligp , ircflp , &
   ischcp , isstpp , inc    , imrgra , ivisse ,                   &
   iwarnp , idftnp , imasac ,                                     &
   blencp , epsrgp , climgp , relaxp , thetap ,                   &
   vel    , vel    ,                                              &
   coefav , coefbv , cofafv , cofbfv ,                            &
   flumas , flumab , viscf  , viscb  , secvif , secvib ,          &
   icvflb , icvfli ,                                              &
   smbr   )

  call field_get_val_s(iestim(iestot), c_estim)
  do iel = 1, ncel
    do isou = 1, 3
      c_estim(iel) = c_estim(iel) + (smbr(isou,iel)/volume(iel))**2
    enddo
  enddo
endif

!===============================================================================
! 4. Finalize the norm of the pressure step (see resopv)
!===============================================================================

if (iappel.eq.1.and.irnpnw.eq.1) then

  ! Compute div(rho u*)

  if (irangp.ge.0.or.iperio.eq.1) then
    call synvin(vel)
  endif

  ! To save time, no space reconstruction
  itypfl = 1
  ! Cavitation algorithm: the pressure step corresponds to the
  ! correction of the volumetric flux, not the mass flux
  if (icavit.ge.0)  itypfl = 0
  init   = 1
  inc    = 1
  iflmb0 = 1
  nswrp  = 1
  imligp = imligr(iu )
  iwarnp = iwarni(ipr)
  epsrgp = epsrgr(iu )
  climgp = climgr(iu )

  call inimav &
 ( ivarfl(iu)      , itypfl ,                                     &
   iflmb0 , init   , inc    , imrgra , nswrp  , imligp ,          &
   iwarnp ,                                                       &
   epsrgp , climgp ,                                              &
   crom   , brom   ,                                              &
   vel    ,                                                       &
   coefav , coefbv ,                                              &
   viscf  , viscb  )

  init = 0
  call divmas(init,viscf,viscb,xnormp)

  ! Compute the norm rnormp used in resopv
  rnormp = sqrt(cs_gdot(ncel,xnormp,xnormp))

  ! Free memory
  deallocate(xnormp)
endif

! ---> Finilaze estimators + Printings

if (iappel.eq.1) then

  ! ---> Estimator on the predicted velocity:
  !      square root (norm) or square root of the sum times the volume (L2 norm)
  if (iestim(iespre).ge.0) then
    call field_get_val_s(iestim(iespre), c_estim)
    if (iescal(iespre).eq.1) then
      do iel = 1, ncel
        c_estim(iel) = sqrt(c_estim(iel))
      enddo
    else if (iescal(iespre).eq.2) then
      do iel = 1, ncel
        c_estim(iel) = sqrt(c_estim(iel)*volume(iel))
      enddo
    endif
  endif

  ! ---> Norm printings
  if (iwarni(iu).ge.2) then
    rnorm = -1.d0
    do iel = 1, ncel
      vitnor = sqrt(vel(1,iel)**2+vel(2,iel)**2+vel(3,iel)**2)
      rnorm = max(rnorm,vitnor)
    enddo

    if (irangp.ge.0) call parmax (rnorm)

    write(nfecra,1100) rnorm

    do iel = 1, ncel
      vitnor = sqrt(vel(1,iel)**2+vel(2,iel)**2+vel(3,iel)**2)
      rnorm = min(rnorm,vitnor)
    enddo

    if (irangp.ge.0) call parmin (rnorm)

    write(nfecra,1200) rnorm

  endif

! ---> Estimator on the whole Navier-Stokes:
!      square root (norm) or square root of the sum times the volume (L2 norm)
else if (iappel.eq.2) then

  call field_get_val_s(iestim(iestot), c_estim)
  if (iescal(iestot).eq.1) then
    do iel = 1, ncel
      c_estim(iel) = sqrt(c_estim(iel))
    enddo
  else if (iescal(iestot).eq.2) then
    do iel = 1, ncel
      c_estim(iel) = sqrt(c_estim(iel)*volume(iel))
    enddo
  endif

endif

! Free memory
!------------
deallocate(smbr)
deallocate(fimp)
deallocate(tsexp)
deallocate(tsimp)
if (allocated(viscce)) deallocate(viscce)
if (allocated(divt)) deallocate(divt)

!--------
! Formats
!--------
#if defined(_CS_LANG_FR)

 1100 format(/,                                                   &
 1X,'Vitesse maximale apres prediction ',E12.4)

 1200 format(/,                                                   &
 1X,'Vitesse minimale apres prediction ',E12.4)

#else

 1100 format(/,                                                   &
 1X,'Maximum velocity after prediction ',E12.4)

 1200 format(/,                                                   &
 1X,'Minimum velocity after prediction ',E12.4)

#endif

!----
! End
!----

return

end subroutine