~ubuntu-branches/ubuntu/precise/code-saturne/precise

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
!-------------------------------------------------------------------------------

! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2011 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.

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

subroutine cptssc &
!================

 ( nvar   , nscal  , ncepdp , ncesmp ,                            &
   iscal  ,                                                       &
   itypfb ,                                                       &
   icepdc , icetsm , itypsm ,                                     &
   izfppp ,                                                       &
   dt     , rtpa   , rtp    , propce , propfa , propfb ,          &
   coefa  , coefb  , ckupdc , smacel ,                            &
   smbrs  , rovsdt )

!===============================================================================
! FONCTION :
! ----------

! ROUTINE PHYSIQUE PARTICULIERE : FLAMME CHARBON PULVERISE
!   ON PRECISE LES TERMES SOURCES POUR UN SCALAIRE PP
!   SUR UN PAS DE TEMPS

! ATTENTION : LE TRAITEMENT DES TERMES SOURCES EST DIFFERENT
! ---------   DE CELUI DE USTSSC.F

! ON RESOUT ROVSDT*D(VAR) = SMBRS

! ROVSDT ET SMBRS CONTIENNENT DEJA D'EVENTUELS TERMES SOURCES
!  UTILISATEUR. IL FAUT DONC LES INCREMENTER ET PAS LES
!  ECRASER

! POUR DES QUESTIONS DE STABILITE, ON NE RAJOUTE DANS ROVSDT
!  QUE DES TERMES POSITIFS. IL N'Y A PAS DE CONTRAINTE POUR
!  SMBRS

! DANS LE CAS D'UN TERME SOURCE EN CEXP + CIMP*VAR ON DOIT
! ECRIRE :
!          SMBRS  = SMBRS  + CEXP + CIMP*VAR
!          ROVSDT = ROVSDT + MAX(-CIMP,ZERO)

! ON FOURNIT ICI ROVSDT ET SMBRS (ILS CONTIENNENT RHO*VOLUME)
!    SMBRS en kg variable/s :
!     ex : pour la vitesse            kg m/s2
!          pour les temperatures      kg degres/s
!          pour les enthalpies        Joules/s
!    ROVSDT en kg /s

!-------------------------------------------------------------------------------
!ARGU                             ARGUMENTS
!__________________.____._____.________________________________________________.
! name             !type!mode ! role                                           !
!__________________!____!_____!________________________________________________!
! nvar             ! i  ! <-- ! total number of variables                      !
! nscal            ! i  ! <-- ! total number of scalars                        !
! ncepdp           ! i  ! <-- ! number of cells with head loss                 !
! ncesmp           ! i  ! <-- ! number of cells with mass source term          !
! iscal            ! i  ! <-- ! scalar number                                  !
! itypfb(nfabor    ! te ! --> ! type des faces de bord                         !
! icepdc(ncelet    ! te ! <-- ! numero des ncepdp cellules avec pdc            !
! icetsm(ncesmp    ! te ! <-- ! numero des cellules a source de masse          !
! itypsm           ! te ! <-- ! type de source de masse pour les               !
! (ncesmp,nvar)    !    !     !  variables (cf. ustsma)                        !
! izfppp           ! te ! --> ! numero de zone de la face de bord              !
! (nfabor)         !    !     !  pour le module phys. part.                    !
! dt(ncelet)       ! ra ! <-- ! time step (per cell)                           !
! rtp, rtpa        ! ra ! <-- ! calculated variables at cell centers           !
!  (ncelet, *)     !    !     !  (at current and previous time steps)          !
! propce(ncelet, *)! ra ! <-- ! physical properties at cell centers            !
! propfa(nfac, *)  ! ra ! <-- ! physical properties at interior face centers   !
! propfb(nfabor, *)! ra ! <-- ! physical properties at boundary face centers   !
! coefa, coefb     ! ra ! <-- ! boundary conditions                            !
!  (nfabor, *)     !    !     !                                                !
! ckupdc           ! tr ! <-- ! tableau de travail pour pdc                    !
!  (ncepdp,6)      !    !     !                                                !
! smacel           ! tr ! <-- ! valeur des variables associee a la             !
! (ncesmp,*   )    !    !     !  source de masse                               !
!                  !    !     !  pour ivar=ipr, smacel=flux de masse           !
! smbrs(ncelet)    ! tr ! --> ! second membre explicite                        !
! rovsdt(ncelet    ! tr ! --> ! partie diagonale implicite                     !
!__________________!____!_____!________________________________________________!

!     TYPE : E (ENTIER), R (REEL), A (ALPHANUMERIQUE), T (TABLEAU)
!            L (LOGIQUE)   .. ET TYPES COMPOSES (EX : TR TABLEAU REEL)
!     MODE : <-- donnee, --> resultat, <-> Donnee modifiee
!            --- tableau de travail
!===============================================================================

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

use paramx
use dimens, only: ndimfb
use numvar
use entsor
use optcal
use cstphy
use cstnum
use parall
use period
use ppppar
use ppthch
use coincl
use cpincl
use ppincl
use ppcpfu
use mesh

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

implicit none

! Arguments

integer          nvar   , nscal
integer          ncepdp , ncesmp
integer          iscal

integer          itypfb(nfabor)
integer          icepdc(ncepdp)
integer          icetsm(ncesmp), itypsm(ncesmp,nvar)
integer          izfppp(nfabor)

double precision dt(ncelet), rtp(ncelet,*), rtpa(ncelet,*)
double precision propce(ncelet,*)
double precision propfa(nfac,*), propfb(ndimfb,*)
double precision coefa(ndimfb,*), coefb(ndimfb,*)
double precision ckupdc(ncepdp,6), smacel(ncesmp,nvar)
double precision smbrs(ncelet), rovsdt(ncelet)

! Local variables

character*80     chaine
integer          ivar , ipcrom, iel
integer          numcla , numcha , icla , numtra
integer          ipcgch , ipcgd1 , ipcgd2 , ipcght , ipcsec
integer          ixchcl , ixckcl , iscala
integer          ipcro2 , ipcte1 , ipcte2 , ipcvsl , ipccp
integer          ipcdia , ipcvst
integer          mode, ige
integer          ifac   , icoefa , icoefb
integer          ipcx2c , icha , imode , ii
integer          iterch
integer          itermx,nbpauv,nbrich,nbepau,nberic,ipghc2

double precision epsrgp , climgp , extrap , xnuss
double precision aux, rhovst
double precision coefe(ngazem)
double precision t1, t2, hlv , hh2ov
double precision f1mc(ncharm), f2mc(ncharm)
double precision xhdev1 , xhdev2 , xhco , xho2 , gamdv1 , gamdv2

double precision gamhet

double precision xxco,xxo2,xxco2,xxh2o,xco2mx
double precision xkp,xk0p,xkm,xk0m,wplus,wmoins,t0p,t0m
double precision auxp,auxm

double precision anmr,xcot,xo2t,xco2e,xo2e,xcoe,tauchi,tautur
double precision sqh2o , x2

double precision err1mx,err2mx,anmr0

double precision errch,ter1,ddelta,xden

double precision, allocatable, dimension(:) :: w1, w2, w3
double precision, allocatable, dimension(:) :: w4, w5, w6
double precision, allocatable, dimension(:) :: w11

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

!===============================================================================
! 1. INITIALISATION
!===============================================================================

! Allocate temporary arrays
allocate(w1(ncelet), w2(ncelet), w3(ncelet))
allocate(w4(ncelet), w5(ncelet), w6(ncelet))
allocate(w11(ncelet))

! Initialize variables to avoid compiler warnings

ipghc2 = 0

! Memoire


! --- Numero du scalaire a traiter : ISCAL

! --- Numero de la variable associee au scalaire a traiter ISCAL
ivar = isca(iscal)

! --- Nom de la variable associee au scalaire a traiter ISCAL
chaine = nomvar(ipprtp(ivar))

! --- Numero des grandeurs physiques (voir usclim)
ipcrom = ipproc(irom)
ipcvst = ipproc(ivisct)

!===============================================================================
! 2. PRISE EN COMPTE DES TERMES SOURCES POUR LES VARIABLES RELATIVES
!    AUX CLASSES DE PARTICULES
!===============================================================================

! --> Terme source pour la fraction massique de charbon reactif

if ( ivar.ge.isca(ixch(1)) .and. ivar.le.isca(ixch(nclacp)) ) then

  if (iwarni(ivar).ge.1) then
    write(nfecra,1000) chaine(1:8)
  endif
  numcla = ivar-isca(ixch(1))+1
  ipcgch = ipproc(igmdch(numcla))

  do iel = 1, ncel

! ---- Calcul de  W1 = - rho.GMDCH > 0

    w1(iel) = - propce(iel,ipcrom)*propce(iel,ipcgch)             &
               *volume(iel)

! ---- Calcul des parties explicite et implicite du TS

    rovsdt(iel) = rovsdt(iel) + max(w1(iel),zero)
    smbrs(iel) = smbrs(iel) - w1(iel)*rtpa(iel,ivar)

  enddo

endif


! --> Terme source pour la fraction massique de coke

if ( ivar.ge.isca(ixck(1)) .and. ivar.le.isca(ixck(nclacp)) ) then

  if (iwarni(ivar).ge.1) then
    write(nfecra,1000) chaine(1:8)
  endif

  numcla = ivar-isca(ixck(1))+1
  ipcgch = ipproc(igmdch(numcla))
  ipcgd1 = ipproc(igmdv1(numcla))
  ipcgd2 = ipproc(igmdv2(numcla))
  ixchcl = isca(ixch(numcla))
  ixckcl = isca(ixck(numcla))
  ipcght = ipproc(igmhet(numcla))
  if ( ihtco2 .eq. 1 ) then
    ipghc2 = ipproc(ighco2(numcla))
  endif

  do iel = 1, ncel

! ---- Calcul de W1 = - rho.Xch.GMDCH.Volume > 0

    w1(iel) = - propce(iel,ipcrom)*rtp(iel,ixchcl)                &
                *propce(iel,ipcgch)*volume(iel)

! AE : On prend RTP(IEL,IXCHCL) et pas RTPA(IEL,IXCHCL) afin
!      d'etre conservatif sur la masse

! ---- Calcul de W2 = rho.Xch.(GMDV1+GMDV2)Volume < 0

    w2(iel) = propce(iel,ipcrom)*rtp(iel,ixchcl)                  &
              *(propce(iel,ipcgd1)+propce(iel,ipcgd2))            &
              *volume(iel)

    if ( rtpa(iel,ixckcl) .gt. epsicp ) then

! Reaction C(s) + O2 ---> 0.5CO
! =============================

! ---- Calcul de la partie implicite  > 0 du TS relatif a GMHET

      w3(iel) = - 2.d0/3.d0*propce(iel,ipcrom)*propce(iel,ipcght) &
               /(rtpa(iel,ixckcl))**(1.d0/3.d0)*volume(iel)

! ---- Calcul de la partie explicite < 0 du TS relatif a GMHET

      w4(iel) = propce(iel,ipcrom)*propce(iel,ipcght)             &
                * (rtpa(iel,ixckcl))**(2.d0/3.d0)*volume(iel)

    else
      w3(iel) = 0.d0
      w4(iel) = 0.d0
    endif

! ---- Calcul des parties explicite et implicite du TS

    rovsdt(iel) = rovsdt(iel) + max(w3(iel),zero)
    smbrs(iel)  = smbrs(iel)  + w1(iel) + w2(iel) + w4(iel)

  enddo

  if ( ihtco2 .eq. 1 ) then

    do iel = 1, ncel

      if ( rtpa(iel,ixckcl) .gt. epsicp ) then

! Reaction C(s) + CO2 ---> 2CO
! =============================

! ---- Calcul de la partie implicite  > 0 du TS relatif a GMHET

        w5(iel) = - 2.d0/3.d0*propce(iel,ipcrom)                  &
                             *propce(iel,ipghc2)                  &
                 /(rtpa(iel,ixckcl))**(1.d0/3.d0)*volume(iel)

! ---- Calcul de la partie explicite < 0 du TS relatif a GMHET

        w6(iel) = propce(iel,ipcrom)*propce(iel,ipghc2)           &
                 *(rtpa(iel,ixckcl))**(2.d0/3.d0)*volume(iel)

      else
        w5(iel) = 0.d0
        w6(iel) = 0.d0
      endif

! ---- Calcul des parties explicite et implicite du TS

      rovsdt(iel) = rovsdt(iel) + max(w5(iel),zero)
      smbrs(iel)  = smbrs(iel)  + w6(iel)

    enddo

  endif

endif


! --> Terme source pour la fraction massique de d'eau

if ( ippmod(icp3pl) .eq. 1 ) then

  if ( ivar.ge.isca(ixwt(1)) .and.                                &
       ivar.le.isca(ixwt(nclacp)) ) then

    if (iwarni(ivar).ge.1) then
      write(nfecra,1000) chaine(1:8)
    endif
    numcla = ivar-isca(ixwt(1))+1
    numcha = ichcor(numcla)

    ipcsec = ipproc(igmsec(numcla))
    ipcx2c = ipproc(ix2(numcla))

    do iel = 1, ncel

! ---- Calcul des parties explicite et implicite du TS

     if ( rtpa(iel,ivar).gt. epsicp .and.                         &
          xwatch(numcha).gt. epsicp       ) then
       w1(iel) = propce(iel,ipcrom)*propce(iel,ipcsec)            &
                     *volume(iel)                                 &
                     *(1.d0/propce(iel,ipcx2c))                   &
                     *(1.d0/xwatch(numcha))

        rovsdt(iel) = rovsdt(iel) + max(w1(iel),zero)
        smbrs(iel)  = smbrs(iel)  - w1(iel)*rtpa(iel,ivar)
     endif

    enddo

  endif

endif

! --> Terme source pour l'enthalpie du solide

if ( ivar.ge.isca(ih2(1)) .and. ivar.le.isca(ih2(nclacp)) ) then

  if (iwarni(ivar).ge.1) then
    write(nfecra,1000) chaine(1:8)
  endif

  numcla = ivar-isca(ih2(1))+1
  numcha = ichcor(numcla)
  ixchcl = isca(ixch(numcla))
  ixckcl = isca(ixck(numcla))
  ipcx2c = ipproc(ix2(numcla))
  ipcro2 = ipproc(irom2(numcla ))
  ipcdia = ipproc(idiam2(numcla))
  ipcte2 = ipproc(itemp2(numcla))
  ipcte1 = ipproc(itemp1)
  ipcght = ipproc(igmhet(numcla))
  if ( ihtco2 .eq. 1 ) then
    ipghc2 = ipproc(ighco2(numcla))
  endif

  ipcgd1 = ipproc(igmdv1(numcla))
  ipcgd2 = ipproc(igmdv2(numcla))
  ipcgch = ipproc(igmdch(numcla))
  ixchcl = isca(ixch(numcla))
  ixckcl = isca(ixck(numcla))

! ---- Calcul preliminaire : calcul de X2(NUMCLA) dans W11
!      ATTENTION tableau a conserver

  do iel = 1, ncel
! Rq : on utilise PROPCE(IEL,IPCX2C) car on veut X2 a l'iteration n
!      (en RTPA)
    w11(iel) = propce(iel,ipcx2c)
  enddo

! ---- Contribution aux bilans explicite et implicite
!        des echanges par diffusion moleculaire
!        6 Lambda Nu / diam**2 / Rho2 * Rho * (T1-T2)

! ------ Calcul de lambda dans W1

  xnuss = 2.d0
  do iel = 1, ncel
    if ( ivisls(ihm).gt.0 ) then
      ipcvsl = ipproc(ivisls(ihm))
      if ( icp.gt.0 ) then
        ipccp   = ipproc(icp)
        w1(iel) = propce(iel,ipcvsl) * propce(iel,ipccp)
      else
        w1(iel) = propce(iel,ipcvsl) * cp0
      endif
    else
      if ( icp.gt.0 ) then
        ipccp   = ipproc(icp)
        w1(iel) = visls0(ihm) * propce(iel,ipccp)
      else
        w1(iel) = visls0(ihm) * cp0
      endif
    endif
  enddo

! ------ Calcul du diametre des particules dans W2
!        On calcule le d20 = (A0.D0**2+(1-A0)*DCK**2)**0.5

  do iel = 1, ncel
    w2(iel) = ( xashch(numcha)*diam20(numcla)**2 +                &
                (1.d0-xashch(numcha))*propce(iel,ipcdia)**2       &
              )**0.5
  enddo

! ------ Contribution aux bilans explicite et implicite de
!        des echanges par diffusion moleculaire

  do iel = 1, ncel
    if ( w11(iel).gt.epsicp ) then
      aux         = 6.d0 * w1(iel) * xnuss / w2(iel)**2           &
                  / propce(iel,ipcro2) * propce(iel,ipcrom)       &
                  * w11(iel)                                      &
                  * volume(iel)
      rhovst      = aux / cp2ch(numcha) /w11(iel)

      smbrs(iel)  = smbrs(iel) - aux *                            &
                  ( propce(iel,ipcte2) - propce(iel,ipcte1) )
      rovsdt(iel) = rovsdt(iel) + max(zero,rhovst)
    endif

  enddo


! ---- Contribution aux bilans explicite et implicite
!        du terme echange d'energie entre les phases :
!        GAMA(dev1) H(mv1,T2)+GAMA(dev2) H(mv2,T2)

  do iel = 1, ncel

!        Gama Dev1 et Gama Dev2

    gamdv1 = propce(iel,ipcrom)*rtp(iel,ixchcl)                   &
            *propce(iel,ipcgd1)

    gamdv2 = propce(iel,ipcrom)*rtp(iel,ixchcl)                   &
            *propce(iel,ipcgd2)

!        H(mv1,T2)

    do ige = 1, ngazem
      coefe(ige) = zero
    enddo

    coefe(ichx1) = a1(numcha)*wmole(ichx1c(numcha))               &
  /( a1(numcha)*wmole(ichx1c(numcha))                             &
    +b1(numcha)*wmole(ico)                                        &
    +c1(numcha)*wmole(ih2o) )
    coefe(ico  ) = b1(numcha)*wmole(ico)                          &
   /( a1(numcha)*wmole(ichx1c(numcha))                            &
     +b1(numcha)*wmole(ico)                                       &
     +c1(numcha)*wmole(ih2o) )
    coefe(ih2o ) = c1(numcha)*wmole(ih2o)                         &
   /( a1(numcha)*wmole(ichx1c(numcha))                            &
     +b1(numcha)*wmole(ico)                                       &
     +c1(numcha)*wmole(ih2o) )

    t2         = propce(iel,ipcte2)
    do icha = 1, ncharm
      f1mc(icha) = zero
      f2mc(icha) = zero
    enddo
    f1mc(numcha) = 1.d0

    mode      = -1
    call cpthp1                                                   &
    !==========
    ( mode  , xhdev1    , coefe  , f1mc   , f2mc   ,              &
      t2    )

!        H(mv2,T2)

    do ige = 1, ngazem
      coefe(ige) = zero
    enddo
    coefe(ichx2) = a2(numcha)*wmole(ichx2c(numcha))               &
   /( a2(numcha)*wmole(ichx2c(numcha))                            &
     +b2(numcha)*wmole(ico)                                       &
     +c2(numcha)*wmole(ih2o) )
    coefe(ico  ) = b2(numcha)*wmole(ico)                          &
   /( a2(numcha)*wmole(ichx2c(numcha))                            &
     +b2(numcha)*wmole(ico)                                       &
     +c2(numcha)*wmole(ih2o) )
    coefe(ih2o ) = c2(numcha)*wmole(ih2o)                         &
   /( a2(numcha)*wmole(ichx2c(numcha))                            &
     +b2(numcha)*wmole(ico)                                       &
     +c2(numcha)*wmole(ih2o) )

    t2         = propce(iel,ipcte2)
    do icha = 1, ncharm
      f1mc(icha) = zero
      f2mc(icha) = zero
    enddo
    f2mc(numcha) = 1.d0

    mode      = -1
    call cpthp1                                                   &
    !==========
    ( mode  , xhdev2    , coefe  , f1mc   , f2mc   ,              &
      t2    )

!         Contribution aux bilans explicite et implicite

    smbrs(iel) = smbrs(iel)                                       &
                 + (gamdv1*xhdev1+gamdv2*xhdev2)*volume(iel)

  enddo

! ------ combustion heterogene : C(s) + 02 ---> 0.5 C0
!        GamHET * (28/12 H(CO,T2)-16/12 H(O2,T1) )

  do iel = 1, ncel

!        Calcul de HCO(T2)

    do ige = 1, ngazem
      coefe(ige) = zero
    enddo
    coefe(ico) = 1.d0
    do icha = 1, ncharm
      f1mc(icha) = zero
      f2mc(icha) = zero
    enddo

    t2        = propce(iel,ipcte2)
    mode      = -1
    call cpthp1                                                   &
    !==========
    ( mode  , xhco    , coefe  , f1mc   , f2mc   ,                &
      t2    )

!        Calcul de HO2(T1)

    do ige = 1, ngazem
      coefe(ige) = zero
    enddo
    coefe(io2) = 1.d0
    do icha = 1, ncharm
      f1mc(icha) = zero
      f2mc(icha) = zero
    enddo

    t1        = propce(iel,ipcte1)
    mode      = -1
    call cpthp1                                                   &
    !==========
    ( mode  , xho2    , coefe  , f1mc   , f2mc   ,                &
      t1    )

!         Contribution aux bilans explicite et implicite

    if ( rtpa(iel,ixckcl) .gt. epsicp ) then

      gamhet = propce(iel,ipcrom)*propce(iel,ipcght)              &
               * ( (rtpa(iel,ixckcl))**(2.d0/3.d0) +              &
                2.d0/3.d0*(rtp(iel,ixckcl)-rtpa(iel,ixckcl))      &
                 /(rtpa(iel,ixckcl))**(1.d0/3.d0) )

    else
      gamhet = 0.d0
    endif

   smbrs(iel) = smbrs(iel)                                        &
                 +gamhet                                          &
                  *(28.d0/12.d0*xhco-16.d0/12.d0*xho2)            &
                  *volume(iel)

  enddo

! ------ combustion heterogene : C(s) + C02 ---> 2 C0
!        GamHET * (56/12 H(CO,T2)-44/12 H(O2,T1) )

  if ( ihtco2 .eq. 1 ) then
    do iel = 1, ncel

!        Calcul de HCO(T2)

      do ige = 1, ngazem
        coefe(ige) = zero
      enddo
      coefe(ico) = 1.d0
      do icha = 1, ncharm
        f1mc(icha) = zero
        f2mc(icha) = zero
      enddo

      t2        = propce(iel,ipcte2)
      mode      = -1
      call cpthp1                                                 &
      !==========
      ( mode  , xhco    , coefe  , f1mc   , f2mc   ,              &
        t2    )

!        Calcul de HCO2(T1)

      do ige = 1, ngazem
        coefe(ige) = zero
      enddo
      coefe(ico2) = 1.d0
      do icha = 1, ncharm
        f1mc(icha) = zero
        f2mc(icha) = zero
      enddo

      t1        = propce(iel,ipcte1)
      mode      = -1
      call cpthp1                                                 &
      !==========
      ( mode  , xho2    , coefe  , f1mc   , f2mc   ,              &
        t1    )

!         Contribution aux bilans explicite et implicite

      if ( rtpa(iel,ixckcl) .gt. epsicp ) then

        gamhet = propce(iel,ipcrom)*propce(iel,ipghc2)            &
                 * ( (rtpa(iel,ixckcl))**(2.d0/3.d0) +            &
                2.d0/3.d0*(rtp(iel,ixckcl)-rtpa(iel,ixckcl))      &
                 /(rtpa(iel,ixckcl))**(1.d0/3.d0) )

      else
        gamhet = 0.d0
      endif

     smbrs(iel) = smbrs(iel)                                      &
                   +gamhet                                        &
                    *(56.d0/12.d0*xhco-44.d0/12.d0*xho2)          &
                    *volume(iel)

    enddo

  endif

!       --> Terme source sur H2 issu du sechage)

  if ( ippmod(icp3pl) .eq. 1 ) then

! ---- Contribution du TS interfacial aux bilans explicite et implicite


    ipcsec = ipproc(igmsec(numcla))
    ipcte2 = ipproc(itemp2(numcla))

    do iel = 1, ncel

!          Calcul de H(H2O) a T2

      do ige = 1, ngazem
        coefe(ige) = zero
      enddo
      coefe(ih2o) = 1.d0
      do icha = 1, ncharm
        f1mc(icha) = zero
        f2mc(icha) = zero
      enddo

      t2 = propce(iel,ipcte2)
      if ( t2 .gt. 100.d0 ) then
        t2 = 100.d0
      endif
      mode      = -1
      call cpthp1                                                 &
      !==========
      ( mode  , hh2ov    , coefe  , f1mc   , f2mc   ,             &
        t2    )

      hlv = 2.263d+6

!         Contribution aux bilans explicite

      if ( rtpa(iel,isca(ixwt(numcla))).gt. epsicp .and.          &
           xwatch(numcha) .gt. epsicp       ) then

        aux = -propce(iel,ipcrom)*propce(iel,ipcsec)              &
       *(rtp(iel,isca(ixwt(numcla)))/propce(iel,ipcx2c))          &
       *(1.d0                    /xwatch(numcha))                 &
             *hh2ov

      else
        aux = 0.d0
      endif

      smbrs(iel) = smbrs(iel) + aux*volume(iel)

    enddo

  endif

endif

!===============================================================================
! 3. PRISE EN COMPTE DES TERMES SOURCES POUR LES VARIABLES RELATIVES
!    AU MELANGE GAZEUX
!===============================================================================

! --> Terme source pour les matieres volatiles legeres

if ( ivar.ge.isca(if1m(1)) .and. ivar.le.isca(if1m(ncharb)) ) then

  if (iwarni(ivar).ge.1) then
    write(nfecra,1000) chaine(1:8)
  endif

! ---- Calcul de GMDEV1 = - SOMME (rho.XCH.GMDV1) > 0  --> W1

  numcha = ivar-isca(if1m(1))+1
  do iel = 1, ncel
    w1(iel) = zero
  enddo
  do icla = 1, nclacp
    ipcgd1 = ipproc(igmdv1(icla))
    ixchcl = isca(ixch(icla))
    if ( ichcor(icla).eq.numcha ) then
      do iel = 1, ncel
        w1(iel) = w1(iel) - propce(iel,ipcrom)*rtp(iel,ixchcl)    &
                * propce(iel,ipcgd1)
      enddo
    endif
  enddo

! ---- Contribution du TS interfacial aux bilans explicite et implicite

  do iel = 1, ncel
    smbrs(iel)  = smbrs(iel)  + volume(iel) * w1(iel)
  enddo

endif


! --> Terme source pour les matieres volatiles lourdes

if ( ivar.ge.isca(if2m(1)) .and. ivar.le.isca(if2m(ncharb)) ) then

  if (iwarni(ivar).ge.1) then
    write(nfecra,1000) chaine(1:8)
  endif

! ---- Calcul de GMDEV2 = - SOMME (rho.XCH.GMDV2) >0 --> W1

  numcha = ivar-isca(if2m(1))+1
  do iel = 1, ncel
    w1(iel) = zero
  enddo
  do icla = 1, nclacp
    ipcgd2 = ipproc(igmdv2(icla))
    ixchcl = isca(ixch(icla))
    if ( ichcor(icla).eq.numcha ) then
      do iel = 1, ncel
        w1(iel) = w1(iel) - propce(iel,ipcrom)*rtp(iel,ixchcl)    &
                * propce(iel,ipcgd2)
      enddo
    endif
  enddo

! ---- Contribution du TS interfacial pour le bilan explicite

  do iel = 1, ncel
    smbrs(iel)  = smbrs(iel)  + volume(iel) * w1(iel)
  enddo

endif


! --> Terme source pour le traceur 3 (O2) (C de la comb. het.)

if ( ivar.eq.isca(if3m) ) then

! RQ IMPORTANTE :  On prend les meme TS que pour Xck
!                  afin d'etre conservatif

  if (iwarni(ivar).ge.1) then
    write(nfecra,1000) chaine(1:8)
  endif

! ---- Calcul prelimimaire pour la partie explicite : W1
!                          pour la partie implicite : W2

  do iel = 1, ncel
    w1(iel) = zero
  enddo

  do icla = 1, nclacp
    ipcght = ipproc(igmhet(icla))
    ixckcl = isca(ixck(icla))
    do iel = 1, ncel
      if ( rtpa(iel,ixckcl) .gt. epsicp ) then
        w1(iel) = w1(iel)                                         &
                 - propce(iel,ipcrom)*propce(iel,ipcght)          &
                 * ( (rtpa(iel,ixckcl))**(2.d0/3.d0) +            &
                    2.d0/3.d0*(rtp(iel,ixckcl)-rtpa(iel,ixckcl))  &
                    /(rtpa(iel,ixckcl))**(1.d0/3.d0) )
      endif
    enddo
  enddo

! ---- Contribution du TS interfacial aux bilans explicite et implicite

  do iel = 1, ncel
    smbrs(iel)  = smbrs(iel)  + volume(iel) * w1(iel)
  enddo

endif



! --> Terme source pour le traceur 3 (CO2) (C de la comb. het.)

if ( ihtco2 .eq. 1 ) then
  if ( ivar.eq.isca(if3mc2) ) then

! RQ IMPORTANTE :  On prend les meme TS que pour Xck
!                  afin d'etre conservatif

    if (iwarni(ivar).ge.1) then
      write(nfecra,1000) chaine(1:8)
    endif

! ---- Calcul prelimimaire pour la partie explicite : W1
!                          pour la partie implicite : W2

    do iel = 1, ncel
      w1(iel) = zero
    enddo

    do icla = 1, nclacp
      ipcght = ipproc(ighco2(icla))
      ixckcl = isca(ixck(icla))
      do iel = 1, ncel
        if ( rtpa(iel,ixckcl) .gt. epsicp ) then
          w1(iel) = w1(iel)                                       &
                   - propce(iel,ipcrom)*propce(iel,ipcght)        &
                   * ( (rtpa(iel,ixckcl))**(2.d0/3.d0) +          &
                      2.d0/3.d0*(rtp(iel,ixckcl)-rtpa(iel,ixckcl))&
                      /(rtpa(iel,ixckcl))**(1.d0/3.d0) )
        endif
      enddo
    enddo

! ---- Contribution du TS interfacial aux bilans explicite et implicite

    do iel = 1, ncel
      smbrs(iel)  = smbrs(iel)  + volume(iel) * w1(iel)
    enddo

  endif

endif

! --> Terme source pour la variance du traceur 4 (Air)

if ( ivar.eq.isca(if4p2m) ) then

  if (iwarni(ivar).ge.1) then
    write(nfecra,1000) chaine(1:8)
  endif

! ---- Calcul des termes sources explicite et implicite
!      relatif aux echanges interfaciaux entre phases

 numtra = 4
 call cptsvi                                                      &
!!==========
 ( ncelet , ncel   , numtra ,                                     &
   rtp    , propce , volume ,                                     &
   smbrs  , rovsdt ,                                              &
   w1     , w2     ,                                              &
   w3 )

! ---- Calcul des termes sources explicite et implicite
!      relatif aux termes de production et de dissipation

!      Pointeur relatif au scalaire associe
!      (0 si pas de scalaire associe)
 iscala = 0

 call cptsvc                                                      &
!!==========
 ( nvar   , nscal  , ncepdp , ncesmp ,                            &
   iscal  , iscala ,          &
   itypfb ,                                                       &
   icepdc , icetsm , itypsm ,                                     &
   dt     , rtpa   , rtp    , propce , propfa , propfb ,          &
   coefa  , coefb  ,                                              &
   smbrs  , rovsdt )

endif

! --> Terme source pour le traceur 5 (Eau issue du séchage)

if ( ippmod(icp3pl) .eq. 1 ) then

  if ( ivar.eq.isca(if5m) ) then


    if (iwarni(ivar).ge.1) then
      write(nfecra,1000) chaine(1:8)
    endif

! ---- Contribution du TS interfacial aux bilans explicite et implicite

    do iel = 1, ncel
      w1(iel) = zero
    enddo

    do icla = 1, nclacp

      ipcsec = ipproc(igmsec(icla))
      ipcx2c = ipproc(ix2(icla))
      numcha = ichcor(icla)

      do iel = 1, ncel

        if (  rtpa(iel,isca(ixwt(icla))).gt. epsicp               &
            .and.                                                 &
              xwatch(numcha) .gt. epsicp       ) then

          w1(iel) = w1(iel)                                       &
      + propce(iel,ipcrom)*propce(iel,ipcsec)                     &
       *(rtp(iel,isca(ixwt(icla)))/propce(iel,ipcx2c))            &
       *(1.d0                  /xwatch(numcha))

        endif

      enddo

    enddo

    do iel = 1, ncel
      smbrs(iel)  = smbrs(iel)  + volume(iel) * w1(iel)
    enddo

  endif

endif

! --> Terme source pour CO2

if ( ieqco2 .ge. 1 ) then

  if ( ivar.eq.isca(iyco2) ) then


    if (iwarni(ivar).ge.1) then
      write(nfecra,1000) chaine(1:8)
    endif

! ---- Contribution du TS interfacial aux bilans explicite et implicite

! Oxydation du CO
! ===============

!  Dryer Glassman : XK0P en (moles/m3)**(-0.75) s-1

!          XK0P = 1.26D10

!           XK0P = 1.26D7 * (1.1)**(NTCABS)
!           IF ( XK0P .GT. 1.26D10 ) XK0P=1.26D10

!           T0P  = 4807.D0

!  Howard : XK0P en (moles/m3)**(-0.75) s-1

!             XK0P = 4.11D9
!             T0P  = 15090.D0

!  Westbrook & Dryer

    xk0p = 23.256d0
    t0p  = 20096.d0

! Dissociation du CO2 (Trinh Minh Chinh)
! ===================

!          XK0M = 5.D8
!          T0M  = 4807.D0

!          XK0M = 0.D0

!  Westbrook & Dryer

    xk0m = 20.03d0
    t0m  = 20096.d0

    err1mx = 0.d0
    err2mx = 0.d0

! Nombre d'iterations

    itermx = 500

! Nombre de points converges

   nbpauv = 0
   nbepau = 0
   nbrich = 0
   nberic = 0

! Precision pour la convergence

   errch = 1.d-8

   do iel = 1, ncel

     xxco  = propce(iel,ipproc(iym1(ico  )))/wmole(ico)           &
            *propce(iel,ipproc(irom1))
     xxo2  = propce(iel,ipproc(iym1(io2  )))/wmole(io2)           &
             *propce(iel,ipproc(irom1))
     xxco2 = propce(iel,ipproc(iym1(ico2 )))/wmole(ico2)          &
            *propce(iel,ipproc(irom1))
     xxh2o = propce(iel,ipproc(iym1(ih2o )))/wmole(ih2o)          &
            *propce(iel,ipproc(irom1))

     xxco  = max(xxco ,zero)
     xxo2  = max(xxo2 ,zero)
     xxco2 = max(xxco2,zero)
     xxh2o = max(xxh2o,zero)

     xco2mx=xxco2+min(xxco,2.d0*xxo2)

!            XKP  = XK0P*EXP(-T0P/PROPCE(IEL,IPPROC(ITEMP1)))
!            XKM  = XK0M*EXP(-T0M/PROPCE(IEL,IPPROC(ITEMP1)))

     xkp = exp(xk0p-t0p/propce(iel,ipproc(itemp1)))
     xkm = exp(xk0m-t0m/propce(iel,ipproc(itemp1)))

!           Recherche de l'état d'équilibre

!           Recerche itérative sans controle de convergence
!            (pour conserver la parallelisation sur les mailles)
!           sur le nombre de moles de reaction séparant
!           l'etat avant réaction (tel que calculé par Cpcym)
!           de l'état d'équilibre

!           initialisation par l'état transporté

     anmr  = xxco2
     xcot  = xxco + xxco2
     xo2t  = xxo2 + 0.5d0*xxco2
     sqh2o = sqrt(xxh2o)

     iterch = 0
 10        continue

       xco2e = anmr
       xcoe  = xcot-anmr
       xo2e  = xo2t-0.5d0*anmr
       iterch = iterch + 1

       ter1 = xkp*(xo2e**0.25d0)*sqh2o
       if ( abs(ter1) .gt. 1.d-15 ) then
         ddelta = (ter1*xcoe-xkm*xco2e)                           &
                 /(ter1*(1.d0+xcoe/(8.d0*xo2e))+xkm)
       else
         ddelta = 0.d0
       endif

       anmr = anmr + ddelta

!  Clipping
!       Cote Pauvre
       if ( xo2t .gt. 0.5d0*xcot ) then
         anmr = max(zero,min(xcot,anmr))
       else
!       Cote Riche
         anmr = max(zero,min(2.d0*xo2t,anmr))
       endif

! Test de convergence

       if ( abs(ddelta) .gt. errch                                &
           .and. iterch.lt.itermx       ) goto 10

! Calcul du nombre de points non convergé

     if ( iterch.eq.itermx                                        &
         .and. abs(ddelta) .gt. errch ) then
       nberic = nberic +1
     endif
     err1mx = max(err1mx,abs(ddelta))

     xco2e = anmr
     xcoe  = xcot-anmr
     xo2e  = xo2t-0.5d0*anmr

     xden = xkp*sqh2o*(0.5d0*(xo2t+xo2e))**0.25d0
     if ( xden .ne. 0.d0 ) then

       tauchi = 1.d0/xden
       tautur = rtpa(iel,ik)/rtpa(iel,iep)

       x2 = 0.d0
       do icla = 1, nclacp
         x2 = x2 + propce(iel,ipproc(ix2(icla)))
       enddo

       if ( ieqco2 .eq. 1 ) then

!              On transporte CO2

         smbrs(iel)  = smbrs(iel)                                 &
                      +wmole(ico2)/propce(iel,ipproc(irom1))      &
         * (xco2e-xxco2)/(tauchi+tautur)                          &
         * (1.d0-x2)                                              &
         * volume(iel) * propce(iel,ipcrom)

       else if ( ieqco2 .eq. 2 ) then

!              On transporte CO

         smbrs(iel)  = smbrs(iel)                                 &
                      +wmole(ico)/propce(iel,ipproc(irom1))       &
         * (xcoe-xxco)/(tauchi+tautur)                            &
         * (1.d0-x2)                                              &
         * volume(iel) * propce(iel,ipcrom)

       endif

       w1(iel) = volume(iel)*propce(iel,ipcrom)/(tauchi+tautur)
       rovsdt(iel) = rovsdt(iel) +   max(w1(iel),zero)

     else
       rovsdt(iel) = rovsdt(iel) + 0.d0
       smbrs(iel)  = smbrs(iel)  + 0.d0
     endif

   enddo

   if(irangp.ge.0) then
     call parcpt(nberic)
     call parmax(err1mx)
   endif
   write(NFECRA,*) ' Erreur max = ',ERR1MX
   write(NFECRA,*) ' Points non   ',NBERIC

!     Terme source : combustion heterogene par le CO2

   if ( ihtco2 .eq. 1) then

     do iel = 1, ncel

       aux = 0.d0
       do icla = 1,nclacp

         ixckcl = isca(ixck(icla))
         ipghc2 = ipproc(ighco2(icla))

         aux = aux                                                &
              + propce(iel,ipcrom)*propce(iel,ipghc2)             &
               *(rtpa(iel,ixckcl))**(2.d0/3.d0)*volume(iel)

       enddo

       rovsdt(iel) = rovsdt(iel) - aux*(wmole(ico2)/0.012)

     enddo

   endif

  endif

endif

! Free memory
deallocate(w1, w2, w3)
deallocate(w4, w5, w6)
deallocate(w11)

!--------
! FORMATS
!--------

 1000 format(' TERMES SOURCES PHYSIQUE PARTICULIERE POUR LA VARIABLE '  &
       ,a8,/)

!----
! FIN
!----

return

end subroutine